14b9ad928SBarry Smith /*$Id: composite.c,v 1.45 2001/08/07 03:03:39 balay Exp $*/ 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines a preconditioner that can consist of a collection of PCs 44b9ad928SBarry Smith */ 54b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 64b9ad928SBarry Smith #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; 124b9ad928SBarry Smith }; 134b9ad928SBarry Smith 144b9ad928SBarry Smith typedef struct { 154b9ad928SBarry Smith PC_CompositeLink head; 164b9ad928SBarry Smith PCCompositeType type; 174b9ad928SBarry Smith Vec work1; 184b9ad928SBarry Smith Vec work2; 194b9ad928SBarry Smith PetscScalar alpha; 204b9ad928SBarry Smith PetscTruth use_true_matrix; 214b9ad928SBarry Smith } PC_Composite; 224b9ad928SBarry Smith 234b9ad928SBarry Smith #undef __FUNCT__ 244b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Multiplicative" 254b9ad928SBarry Smith static int PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y) 264b9ad928SBarry Smith { 274b9ad928SBarry Smith int ierr; 284b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 294b9ad928SBarry Smith PC_CompositeLink next = jac->head; 304b9ad928SBarry Smith PetscScalar one = 1.0,mone = -1.0; 314b9ad928SBarry Smith Mat mat = pc->pmat; 324b9ad928SBarry Smith 334b9ad928SBarry Smith PetscFunctionBegin; 344b9ad928SBarry Smith if (!next) { 354b9ad928SBarry Smith SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 364b9ad928SBarry Smith } 374b9ad928SBarry Smith if (next->next && !jac->work2) { /* allocate second work vector */ 384b9ad928SBarry Smith ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr); 394b9ad928SBarry Smith } 40*d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 414b9ad928SBarry Smith if (jac->use_true_matrix) mat = pc->mat; 424b9ad928SBarry Smith while (next->next) { 434b9ad928SBarry Smith next = next->next; 444b9ad928SBarry Smith ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr); 454b9ad928SBarry Smith ierr = VecWAXPY(&mone,jac->work1,x,jac->work2);CHKERRQ(ierr); 46*d8fd42c4SBarry Smith ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 474b9ad928SBarry Smith ierr = VecAXPY(&one,jac->work1,y);CHKERRQ(ierr); 484b9ad928SBarry Smith } 494b9ad928SBarry Smith 504b9ad928SBarry Smith PetscFunctionReturn(0); 514b9ad928SBarry Smith } 524b9ad928SBarry Smith 534b9ad928SBarry Smith /* 544b9ad928SBarry Smith This is very special for a matrix of the form alpha I + R + S 554b9ad928SBarry Smith where first preconditioner is built from alpha I + S and second from 564b9ad928SBarry Smith alpha I + R 574b9ad928SBarry Smith */ 584b9ad928SBarry Smith #undef __FUNCT__ 594b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Special" 604b9ad928SBarry Smith static int PCApply_Composite_Special(PC pc,Vec x,Vec y) 614b9ad928SBarry Smith { 624b9ad928SBarry Smith int ierr; 634b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 644b9ad928SBarry Smith PC_CompositeLink next = jac->head; 654b9ad928SBarry Smith 664b9ad928SBarry Smith PetscFunctionBegin; 674b9ad928SBarry Smith if (!next) { 684b9ad928SBarry Smith SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 694b9ad928SBarry Smith } 704b9ad928SBarry Smith if (!next->next || next->next->next) { 714b9ad928SBarry Smith SETERRQ(1,"Special composite preconditioners requires exactly two PCs"); 724b9ad928SBarry Smith } 734b9ad928SBarry Smith 74*d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 75*d8fd42c4SBarry Smith ierr = PCApply(next->next->pc,jac->work1,y);CHKERRQ(ierr); 764b9ad928SBarry Smith PetscFunctionReturn(0); 774b9ad928SBarry Smith } 784b9ad928SBarry Smith 794b9ad928SBarry Smith #undef __FUNCT__ 804b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Additive" 814b9ad928SBarry Smith static int PCApply_Composite_Additive(PC pc,Vec x,Vec y) 824b9ad928SBarry Smith { 834b9ad928SBarry Smith int ierr; 844b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 854b9ad928SBarry Smith PC_CompositeLink next = jac->head; 864b9ad928SBarry Smith PetscScalar one = 1.0; 874b9ad928SBarry Smith 884b9ad928SBarry Smith PetscFunctionBegin; 894b9ad928SBarry Smith if (!next) { 904b9ad928SBarry Smith SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 914b9ad928SBarry Smith } 92*d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 934b9ad928SBarry Smith while (next->next) { 944b9ad928SBarry Smith next = next->next; 95*d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 964b9ad928SBarry Smith ierr = VecAXPY(&one,jac->work1,y);CHKERRQ(ierr); 974b9ad928SBarry Smith } 984b9ad928SBarry Smith PetscFunctionReturn(0); 994b9ad928SBarry Smith } 1004b9ad928SBarry Smith 1014b9ad928SBarry Smith #undef __FUNCT__ 1024b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Composite" 1034b9ad928SBarry Smith static int PCSetUp_Composite(PC pc) 1044b9ad928SBarry Smith { 1054b9ad928SBarry Smith int ierr; 1064b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1074b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1084b9ad928SBarry Smith 1094b9ad928SBarry Smith PetscFunctionBegin; 1104b9ad928SBarry Smith if (!jac->work1) { 11123ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr); 1124b9ad928SBarry Smith } 1134b9ad928SBarry Smith while (next) { 1144b9ad928SBarry Smith ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1154b9ad928SBarry Smith next = next->next; 1164b9ad928SBarry Smith } 1174b9ad928SBarry Smith 1184b9ad928SBarry Smith PetscFunctionReturn(0); 1194b9ad928SBarry Smith } 1204b9ad928SBarry Smith 1214b9ad928SBarry Smith #undef __FUNCT__ 1224b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Composite" 1234b9ad928SBarry Smith static int PCDestroy_Composite(PC pc) 1244b9ad928SBarry Smith { 1254b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1264b9ad928SBarry Smith int ierr; 1274b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith PetscFunctionBegin; 1304b9ad928SBarry Smith while (next) { 1314b9ad928SBarry Smith ierr = PCDestroy(next->pc);CHKERRQ(ierr); 1324b9ad928SBarry Smith next = next->next; 1334b9ad928SBarry Smith } 1344b9ad928SBarry Smith 1354b9ad928SBarry Smith if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);} 1364b9ad928SBarry Smith if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);} 1374b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 1384b9ad928SBarry Smith PetscFunctionReturn(0); 1394b9ad928SBarry Smith } 1404b9ad928SBarry Smith 1414b9ad928SBarry Smith #undef __FUNCT__ 1424b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite" 1434b9ad928SBarry Smith static int PCSetFromOptions_Composite(PC pc) 1444b9ad928SBarry Smith { 1454b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1464b9ad928SBarry Smith int ierr,nmax = 8,i,indx; 1474b9ad928SBarry Smith PC_CompositeLink next; 148e5999256SBarry Smith char *pcs[8]; 149e5999256SBarry Smith const char *types[] = {"multiplicative","additive","special"}; 1504b9ad928SBarry Smith PetscTruth flg; 1514b9ad928SBarry Smith 1524b9ad928SBarry Smith PetscFunctionBegin; 1534b9ad928SBarry Smith ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 1544b9ad928SBarry Smith ierr = PetscOptionsEList("-pc_composite_type","Type of composition","PCCompositeSetType",types,3,types[0],&indx,&flg);CHKERRQ(ierr); 1554b9ad928SBarry Smith if (flg) { 1564b9ad928SBarry Smith ierr = PCCompositeSetType(pc,(PCCompositeType)indx);CHKERRQ(ierr); 1574b9ad928SBarry Smith } 1584b9ad928SBarry Smith ierr = PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);CHKERRQ(ierr); 1594b9ad928SBarry Smith if (flg) { 1604b9ad928SBarry Smith ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr); 1614b9ad928SBarry Smith } 1624b9ad928SBarry Smith ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 1634b9ad928SBarry Smith if (flg) { 1644b9ad928SBarry Smith for (i=0; i<nmax; i++) { 1654b9ad928SBarry Smith ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 1664b9ad928SBarry Smith } 1674b9ad928SBarry Smith } 1684b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 1694b9ad928SBarry Smith 1704b9ad928SBarry Smith next = jac->head; 1714b9ad928SBarry Smith while (next) { 1724b9ad928SBarry Smith ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 1734b9ad928SBarry Smith next = next->next; 1744b9ad928SBarry Smith } 1754b9ad928SBarry Smith PetscFunctionReturn(0); 1764b9ad928SBarry Smith } 1774b9ad928SBarry Smith 1784b9ad928SBarry Smith #undef __FUNCT__ 1794b9ad928SBarry Smith #define __FUNCT__ "PCView_Composite" 1804b9ad928SBarry Smith static int PCView_Composite(PC pc,PetscViewer viewer) 1814b9ad928SBarry Smith { 1824b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1834b9ad928SBarry Smith int ierr; 1844b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1854b9ad928SBarry Smith PetscTruth isascii; 1864b9ad928SBarry Smith 1874b9ad928SBarry Smith PetscFunctionBegin; 1884b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 1894b9ad928SBarry Smith if (isascii) { 1904b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 1914b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 1924b9ad928SBarry Smith } else { 1934b9ad928SBarry Smith SETERRQ1(1,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name); 1944b9ad928SBarry Smith } 1954b9ad928SBarry Smith if (isascii) { 1964b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1974b9ad928SBarry Smith } 1984b9ad928SBarry Smith while (next) { 1994b9ad928SBarry Smith ierr = PCView(next->pc,viewer);CHKERRQ(ierr); 2004b9ad928SBarry Smith next = next->next; 2014b9ad928SBarry Smith } 2024b9ad928SBarry Smith if (isascii) { 2034b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2044b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 2054b9ad928SBarry Smith } 2064b9ad928SBarry Smith PetscFunctionReturn(0); 2074b9ad928SBarry Smith } 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 2104b9ad928SBarry Smith 2114b9ad928SBarry Smith EXTERN_C_BEGIN 2124b9ad928SBarry Smith #undef __FUNCT__ 2134b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite" 2144b9ad928SBarry Smith int PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 2154b9ad928SBarry Smith { 2164b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2174b9ad928SBarry Smith PetscFunctionBegin; 2184b9ad928SBarry Smith jac->alpha = alpha; 2194b9ad928SBarry Smith PetscFunctionReturn(0); 2204b9ad928SBarry Smith } 2214b9ad928SBarry Smith EXTERN_C_END 2224b9ad928SBarry Smith 2234b9ad928SBarry Smith EXTERN_C_BEGIN 2244b9ad928SBarry Smith #undef __FUNCT__ 2254b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType_Composite" 2264b9ad928SBarry Smith int PCCompositeSetType_Composite(PC pc,PCCompositeType type) 2274b9ad928SBarry Smith { 2284b9ad928SBarry Smith PetscFunctionBegin; 2294b9ad928SBarry Smith if (type == PC_COMPOSITE_ADDITIVE) { 2304b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 2314b9ad928SBarry Smith } else if (type == PC_COMPOSITE_MULTIPLICATIVE) { 2324b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Multiplicative; 2334b9ad928SBarry Smith } else if (type == PC_COMPOSITE_SPECIAL) { 2344b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Special; 2354b9ad928SBarry Smith } else { 2364b9ad928SBarry Smith SETERRQ(1,"Unkown composite preconditioner type"); 2374b9ad928SBarry Smith } 2384b9ad928SBarry Smith PetscFunctionReturn(0); 2394b9ad928SBarry Smith } 2404b9ad928SBarry Smith EXTERN_C_END 2414b9ad928SBarry Smith 2424b9ad928SBarry Smith EXTERN_C_BEGIN 2434b9ad928SBarry Smith #undef __FUNCT__ 2444b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC_Composite" 2454b9ad928SBarry Smith int PCCompositeAddPC_Composite(PC pc,PCType type) 2464b9ad928SBarry Smith { 2474b9ad928SBarry Smith PC_Composite *jac; 2484b9ad928SBarry Smith PC_CompositeLink next,link; 2494b9ad928SBarry Smith int ierr,cnt = 0; 2504b9ad928SBarry Smith char *prefix,newprefix[8]; 2514b9ad928SBarry Smith 2524b9ad928SBarry Smith PetscFunctionBegin; 2534b9ad928SBarry Smith ierr = PetscNew(struct _PC_CompositeLink,&link);CHKERRQ(ierr); 2544b9ad928SBarry Smith link->next = 0; 2554b9ad928SBarry Smith ierr = PCCreate(pc->comm,&link->pc);CHKERRQ(ierr); 2564b9ad928SBarry Smith 2574b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 2584b9ad928SBarry Smith next = jac->head; 2594b9ad928SBarry Smith if (!next) { 2604b9ad928SBarry Smith jac->head = link; 2614b9ad928SBarry Smith } else { 2624b9ad928SBarry Smith cnt++; 2634b9ad928SBarry Smith while (next->next) { 2644b9ad928SBarry Smith next = next->next; 2654b9ad928SBarry Smith cnt++; 2664b9ad928SBarry Smith } 2674b9ad928SBarry Smith next->next = link; 2684b9ad928SBarry Smith } 2694b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 2704b9ad928SBarry Smith ierr = PCSetOptionsPrefix(link->pc,prefix);CHKERRQ(ierr); 2714b9ad928SBarry Smith sprintf(newprefix,"sub_%d_",cnt); 2724b9ad928SBarry Smith ierr = PCAppendOptionsPrefix(link->pc,newprefix);CHKERRQ(ierr); 2734b9ad928SBarry Smith /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 2744b9ad928SBarry Smith ierr = PCSetType(link->pc,type);CHKERRQ(ierr); 2754b9ad928SBarry Smith 2764b9ad928SBarry Smith PetscFunctionReturn(0); 2774b9ad928SBarry Smith } 2784b9ad928SBarry Smith EXTERN_C_END 2794b9ad928SBarry Smith 2804b9ad928SBarry Smith EXTERN_C_BEGIN 2814b9ad928SBarry Smith #undef __FUNCT__ 2824b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC_Composite" 2834b9ad928SBarry Smith int PCCompositeGetPC_Composite(PC pc,int n,PC *subpc) 2844b9ad928SBarry Smith { 2854b9ad928SBarry Smith PC_Composite *jac; 2864b9ad928SBarry Smith PC_CompositeLink next; 2874b9ad928SBarry Smith int i; 2884b9ad928SBarry Smith 2894b9ad928SBarry Smith PetscFunctionBegin; 2904b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 2914b9ad928SBarry Smith next = jac->head; 2924b9ad928SBarry Smith for (i=0; i<n; i++) { 2934b9ad928SBarry Smith if (!next->next) { 2944b9ad928SBarry Smith SETERRQ(1,"Not enough PCs in composite preconditioner"); 2954b9ad928SBarry Smith } 2964b9ad928SBarry Smith next = next->next; 2974b9ad928SBarry Smith } 2984b9ad928SBarry Smith *subpc = next->pc; 2994b9ad928SBarry Smith PetscFunctionReturn(0); 3004b9ad928SBarry Smith } 3014b9ad928SBarry Smith EXTERN_C_END 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith EXTERN_C_BEGIN 3044b9ad928SBarry Smith #undef __FUNCT__ 3054b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue_Composite" 3064b9ad928SBarry Smith int PCCompositeSetUseTrue_Composite(PC pc) 3074b9ad928SBarry Smith { 3084b9ad928SBarry Smith PC_Composite *jac; 3094b9ad928SBarry Smith 3104b9ad928SBarry Smith PetscFunctionBegin; 3114b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 3124b9ad928SBarry Smith jac->use_true_matrix = PETSC_TRUE; 3134b9ad928SBarry Smith PetscFunctionReturn(0); 3144b9ad928SBarry Smith } 3154b9ad928SBarry Smith EXTERN_C_END 3164b9ad928SBarry Smith 3174b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */ 3184b9ad928SBarry Smith #undef __FUNCT__ 3194b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType" 3204b9ad928SBarry Smith /*@C 3214b9ad928SBarry Smith PCCompositeSetType - Sets the type of composite preconditioner. 3224b9ad928SBarry Smith 3234b9ad928SBarry Smith Collective on PC 3244b9ad928SBarry Smith 3254b9ad928SBarry Smith Input Parameter: 3264b9ad928SBarry Smith . pc - the preconditioner context 3274b9ad928SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 3284b9ad928SBarry Smith 3294b9ad928SBarry Smith Options Database Key: 3304b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith Level: Developer 3334b9ad928SBarry Smith 3344b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 3354b9ad928SBarry Smith @*/ 3364b9ad928SBarry Smith int PCCompositeSetType(PC pc,PCCompositeType type) 3374b9ad928SBarry Smith { 3384b9ad928SBarry Smith int ierr,(*f)(PC,PCCompositeType); 3394b9ad928SBarry Smith 3404b9ad928SBarry Smith PetscFunctionBegin; 3414482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3424b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 3434b9ad928SBarry Smith if (f) { 3444b9ad928SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 3454b9ad928SBarry Smith } 3464b9ad928SBarry Smith PetscFunctionReturn(0); 3474b9ad928SBarry Smith } 3484b9ad928SBarry Smith 3494b9ad928SBarry Smith #undef __FUNCT__ 3504b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha" 3514b9ad928SBarry Smith /*@C 3524b9ad928SBarry Smith PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner 3534b9ad928SBarry Smith for alphaI + R + S 3544b9ad928SBarry Smith 3554b9ad928SBarry Smith Collective on PC 3564b9ad928SBarry Smith 3574b9ad928SBarry Smith Input Parameter: 3584b9ad928SBarry Smith + pc - the preconditioner context 3594b9ad928SBarry Smith - alpha - scale on identity 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith Level: Developer 3624b9ad928SBarry Smith 3634b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 3644b9ad928SBarry Smith @*/ 3654b9ad928SBarry Smith int PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 3664b9ad928SBarry Smith { 3674b9ad928SBarry Smith int ierr,(*f)(PC,PetscScalar); 3684b9ad928SBarry Smith 3694b9ad928SBarry Smith PetscFunctionBegin; 3704482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3714b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);CHKERRQ(ierr); 3724b9ad928SBarry Smith if (f) { 3734b9ad928SBarry Smith ierr = (*f)(pc,alpha);CHKERRQ(ierr); 3744b9ad928SBarry Smith } 3754b9ad928SBarry Smith PetscFunctionReturn(0); 3764b9ad928SBarry Smith } 3774b9ad928SBarry Smith 3784b9ad928SBarry Smith #undef __FUNCT__ 3794b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC" 3804b9ad928SBarry Smith /*@C 3814b9ad928SBarry Smith PCCompositeAddPC - Adds another PC to the composite PC. 3824b9ad928SBarry Smith 3834b9ad928SBarry Smith Collective on PC 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith Input Parameters: 3864b9ad928SBarry Smith . pc - the preconditioner context 3874b9ad928SBarry Smith . type - the type of the new preconditioner 3884b9ad928SBarry Smith 3894b9ad928SBarry Smith Level: Developer 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith .keywords: PC, composite preconditioner, add 3924b9ad928SBarry Smith @*/ 3934b9ad928SBarry Smith int PCCompositeAddPC(PC pc,PCType type) 3944b9ad928SBarry Smith { 3954b9ad928SBarry Smith int ierr,(*f)(PC,PCType); 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith PetscFunctionBegin; 3984482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3994b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);CHKERRQ(ierr); 4004b9ad928SBarry Smith if (f) { 4014b9ad928SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 4024b9ad928SBarry Smith } 4034b9ad928SBarry Smith PetscFunctionReturn(0); 4044b9ad928SBarry Smith } 4054b9ad928SBarry Smith 4064b9ad928SBarry Smith #undef __FUNCT__ 4074b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC" 4084b9ad928SBarry Smith /*@C 4094b9ad928SBarry Smith PCCompositeGetPC - Gets one of the PC objects in the composite PC. 4104b9ad928SBarry Smith 4114b9ad928SBarry Smith Not Collective 4124b9ad928SBarry Smith 4134b9ad928SBarry Smith Input Parameter: 4144b9ad928SBarry Smith . pc - the preconditioner context 4154b9ad928SBarry Smith . n - the number of the pc requested 4164b9ad928SBarry Smith 4174b9ad928SBarry Smith Output Parameters: 4184b9ad928SBarry Smith . subpc - the PC requested 4194b9ad928SBarry Smith 4204b9ad928SBarry Smith Level: Developer 4214b9ad928SBarry Smith 4224b9ad928SBarry Smith .keywords: PC, get, composite preconditioner, sub preconditioner 4234b9ad928SBarry Smith 4244b9ad928SBarry Smith .seealso: PCCompositeAddPC() 4254b9ad928SBarry Smith @*/ 4264b9ad928SBarry Smith int PCCompositeGetPC(PC pc,int n,PC *subpc) 4274b9ad928SBarry Smith { 4284b9ad928SBarry Smith int ierr,(*f)(PC,int,PC *); 4294b9ad928SBarry Smith 4304b9ad928SBarry Smith PetscFunctionBegin; 4314482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4324482741eSBarry Smith PetscValidPointer(subpc,3); 4334b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 4344b9ad928SBarry Smith if (f) { 4354b9ad928SBarry Smith ierr = (*f)(pc,n,subpc);CHKERRQ(ierr); 4364b9ad928SBarry Smith } else { 4374b9ad928SBarry Smith SETERRQ(1,"Cannot get pc, not composite type"); 4384b9ad928SBarry Smith } 4394b9ad928SBarry Smith PetscFunctionReturn(0); 4404b9ad928SBarry Smith } 4414b9ad928SBarry Smith 4424b9ad928SBarry Smith #undef __FUNCT__ 4434b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue" 4444b9ad928SBarry Smith /*@ 4454b9ad928SBarry Smith PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than 4464b9ad928SBarry Smith the matrix used to define the preconditioner) is used to compute 4474b9ad928SBarry Smith the residual when the multiplicative scheme is used. 4484b9ad928SBarry Smith 4494b9ad928SBarry Smith Collective on PC 4504b9ad928SBarry Smith 4514b9ad928SBarry Smith Input Parameters: 4524b9ad928SBarry Smith . pc - the preconditioner context 4534b9ad928SBarry Smith 4544b9ad928SBarry Smith Options Database Key: 4554b9ad928SBarry Smith . -pc_composite_true - Activates PCCompositeSetUseTrue() 4564b9ad928SBarry Smith 4574b9ad928SBarry Smith Note: 4584b9ad928SBarry Smith For the common case in which the preconditioning and linear 4594b9ad928SBarry Smith system matrices are identical, this routine is unnecessary. 4604b9ad928SBarry Smith 4614b9ad928SBarry Smith Level: Developer 4624b9ad928SBarry Smith 4634b9ad928SBarry Smith .keywords: PC, composite preconditioner, set, true, flag 4644b9ad928SBarry Smith 4654b9ad928SBarry Smith .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue() 4664b9ad928SBarry Smith @*/ 4674b9ad928SBarry Smith int PCCompositeSetUseTrue(PC pc) 4684b9ad928SBarry Smith { 4694b9ad928SBarry Smith int ierr,(*f)(PC); 4704b9ad928SBarry Smith 4714b9ad928SBarry Smith PetscFunctionBegin; 4724482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4734b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr); 4744b9ad928SBarry Smith if (f) { 4754b9ad928SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 4764b9ad928SBarry Smith } 4774b9ad928SBarry Smith PetscFunctionReturn(0); 4784b9ad928SBarry Smith } 4794b9ad928SBarry Smith 4804b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/ 4814b9ad928SBarry Smith 4824b9ad928SBarry Smith /*MC 4834b9ad928SBarry Smith PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 4844b9ad928SBarry Smith 4854b9ad928SBarry Smith Options Database Keys: 4864b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 4874b9ad928SBarry Smith . -pc_composite_true - Activates PCCompositeSetUseTrue() 4884b9ad928SBarry Smith 4894b9ad928SBarry Smith Level: intermediate 4904b9ad928SBarry Smith 4914b9ad928SBarry Smith Concepts: composing solvers 4924b9ad928SBarry Smith 4934b9ad928SBarry Smith Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 4944b9ad928SBarry Smith inner PCs to be PCKSP. 4954b9ad928SBarry Smith Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 4964b9ad928SBarry Smith the incorrect answer) unless you use KSPFGMRES as the other Krylov method 4974b9ad928SBarry Smith 4984b9ad928SBarry Smith 4994b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 5004b9ad928SBarry Smith PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 5014b9ad928SBarry Smith PCCompositeGetPC(), PCCompositeSetUseTrue() 5024b9ad928SBarry Smith 5034b9ad928SBarry Smith M*/ 5044b9ad928SBarry Smith 5054b9ad928SBarry Smith EXTERN_C_BEGIN 5064b9ad928SBarry Smith #undef __FUNCT__ 5074b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Composite" 5084b9ad928SBarry Smith int PCCreate_Composite(PC pc) 5094b9ad928SBarry Smith { 5104b9ad928SBarry Smith int ierr; 5114b9ad928SBarry Smith PC_Composite *jac; 5124b9ad928SBarry Smith 5134b9ad928SBarry Smith PetscFunctionBegin; 5144b9ad928SBarry Smith ierr = PetscNew(PC_Composite,&jac);CHKERRQ(ierr); 5154b9ad928SBarry Smith PetscLogObjectMemory(pc,sizeof(PC_Composite)); 5164b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 5174b9ad928SBarry Smith pc->ops->setup = PCSetUp_Composite; 5184b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Composite; 5194b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Composite; 5204b9ad928SBarry Smith pc->ops->view = PCView_Composite; 5214b9ad928SBarry Smith pc->ops->applyrichardson = 0; 5224b9ad928SBarry Smith 5234b9ad928SBarry Smith pc->data = (void*)jac; 5244b9ad928SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 5254b9ad928SBarry Smith jac->work1 = 0; 5264b9ad928SBarry Smith jac->work2 = 0; 5274b9ad928SBarry Smith jac->head = 0; 5284b9ad928SBarry Smith 5294b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite", 5304b9ad928SBarry Smith PCCompositeSetType_Composite);CHKERRQ(ierr); 5314b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite", 5324b9ad928SBarry Smith PCCompositeAddPC_Composite);CHKERRQ(ierr); 5334b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite", 5344b9ad928SBarry Smith PCCompositeGetPC_Composite);CHKERRQ(ierr); 5354b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite", 5364b9ad928SBarry Smith PCCompositeSetUseTrue_Composite);CHKERRQ(ierr); 5374b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite", 5384b9ad928SBarry Smith PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr); 5394b9ad928SBarry Smith 5404b9ad928SBarry Smith PetscFunctionReturn(0); 5414b9ad928SBarry Smith } 5424b9ad928SBarry Smith EXTERN_C_END 5434b9ad928SBarry Smith 544