1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 34b9ad928SBarry Smith /* 44b9ad928SBarry Smith Defines a preconditioner that can consist of a collection of PCs 54b9ad928SBarry Smith */ 64b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 74b9ad928SBarry Smith #include "petscksp.h" /*I "petscksp.h" I*/ 84b9ad928SBarry Smith 94b9ad928SBarry Smith typedef struct _PC_CompositeLink *PC_CompositeLink; 104b9ad928SBarry Smith struct _PC_CompositeLink { 114b9ad928SBarry Smith PC pc; 124b9ad928SBarry Smith PC_CompositeLink next; 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 PetscTruth use_true_matrix; 224b9ad928SBarry Smith } PC_Composite; 234b9ad928SBarry Smith 244b9ad928SBarry Smith #undef __FUNCT__ 254b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Multiplicative" 266849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y) 274b9ad928SBarry Smith { 28dfbe8321SBarry Smith PetscErrorCode ierr; 294b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 304b9ad928SBarry Smith PC_CompositeLink next = jac->head; 314b9ad928SBarry Smith PetscScalar one = 1.0,mone = -1.0; 324b9ad928SBarry Smith Mat mat = pc->pmat; 334b9ad928SBarry Smith 344b9ad928SBarry Smith PetscFunctionBegin; 354b9ad928SBarry Smith if (!next) { 361302d50aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()"); 374b9ad928SBarry Smith } 384b9ad928SBarry Smith if (next->next && !jac->work2) { /* allocate second work vector */ 394b9ad928SBarry Smith ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr); 404b9ad928SBarry Smith } 41d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 424b9ad928SBarry Smith if (jac->use_true_matrix) mat = pc->mat; 434b9ad928SBarry Smith while (next->next) { 444b9ad928SBarry Smith next = next->next; 454b9ad928SBarry Smith ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr); 462dcb1b2aSMatthew Knepley ierr = VecWAXPY(jac->work2,mone,jac->work1,x);CHKERRQ(ierr); 47d8fd42c4SBarry Smith ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 482dcb1b2aSMatthew Knepley ierr = VecAXPY(y,one,jac->work1);CHKERRQ(ierr); 494b9ad928SBarry Smith } 504b9ad928SBarry Smith 514b9ad928SBarry Smith PetscFunctionReturn(0); 524b9ad928SBarry Smith } 534b9ad928SBarry Smith 544b9ad928SBarry Smith /* 554b9ad928SBarry Smith This is very special for a matrix of the form alpha I + R + S 564b9ad928SBarry Smith where first preconditioner is built from alpha I + S and second from 574b9ad928SBarry Smith alpha I + R 584b9ad928SBarry Smith */ 594b9ad928SBarry Smith #undef __FUNCT__ 604b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Special" 616849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y) 624b9ad928SBarry Smith { 63dfbe8321SBarry Smith PetscErrorCode ierr; 644b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 654b9ad928SBarry Smith PC_CompositeLink next = jac->head; 664b9ad928SBarry Smith 674b9ad928SBarry Smith PetscFunctionBegin; 684b9ad928SBarry Smith if (!next) { 691302d50aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()"); 704b9ad928SBarry Smith } 714b9ad928SBarry Smith if (!next->next || next->next->next) { 721302d50aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Special composite preconditioners requires exactly two PCs"); 734b9ad928SBarry Smith } 744b9ad928SBarry Smith 75d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 76d8fd42c4SBarry Smith ierr = PCApply(next->next->pc,jac->work1,y);CHKERRQ(ierr); 774b9ad928SBarry Smith PetscFunctionReturn(0); 784b9ad928SBarry Smith } 794b9ad928SBarry Smith 804b9ad928SBarry Smith #undef __FUNCT__ 814b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Additive" 826849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y) 834b9ad928SBarry Smith { 84dfbe8321SBarry Smith PetscErrorCode ierr; 854b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 864b9ad928SBarry Smith PC_CompositeLink next = jac->head; 874b9ad928SBarry Smith PetscScalar one = 1.0; 884b9ad928SBarry Smith 894b9ad928SBarry Smith PetscFunctionBegin; 904b9ad928SBarry Smith if (!next) { 911302d50aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No composite preconditioners supplied via PCCompositeAddPC()"); 924b9ad928SBarry Smith } 93d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 944b9ad928SBarry Smith while (next->next) { 954b9ad928SBarry Smith next = next->next; 96d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 972dcb1b2aSMatthew Knepley ierr = VecAXPY(y,one,jac->work1);CHKERRQ(ierr); 984b9ad928SBarry Smith } 994b9ad928SBarry Smith PetscFunctionReturn(0); 1004b9ad928SBarry Smith } 1014b9ad928SBarry Smith 1024b9ad928SBarry Smith #undef __FUNCT__ 1034b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Composite" 1046849ba73SBarry Smith static PetscErrorCode PCSetUp_Composite(PC pc) 1054b9ad928SBarry Smith { 106dfbe8321SBarry Smith PetscErrorCode ierr; 1074b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1084b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1094b9ad928SBarry Smith 1104b9ad928SBarry Smith PetscFunctionBegin; 1114b9ad928SBarry Smith if (!jac->work1) { 11223ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr); 1134b9ad928SBarry Smith } 1144b9ad928SBarry Smith while (next) { 1154b9ad928SBarry Smith ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1164b9ad928SBarry Smith next = next->next; 1174b9ad928SBarry Smith } 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith PetscFunctionReturn(0); 1204b9ad928SBarry Smith } 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith #undef __FUNCT__ 1234b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Composite" 1246849ba73SBarry Smith static PetscErrorCode PCDestroy_Composite(PC pc) 1254b9ad928SBarry Smith { 1264b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 127dfbe8321SBarry Smith PetscErrorCode ierr; 1284b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1294b9ad928SBarry Smith 1304b9ad928SBarry Smith PetscFunctionBegin; 1314b9ad928SBarry Smith while (next) { 1324b9ad928SBarry Smith ierr = PCDestroy(next->pc);CHKERRQ(ierr); 1334b9ad928SBarry Smith next = next->next; 1344b9ad928SBarry Smith } 1354b9ad928SBarry Smith 1364b9ad928SBarry Smith if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);} 1374b9ad928SBarry Smith if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);} 1384b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 1394b9ad928SBarry Smith PetscFunctionReturn(0); 1404b9ad928SBarry Smith } 1414b9ad928SBarry Smith 1424b9ad928SBarry Smith #undef __FUNCT__ 1434b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite" 1446849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PC pc) 1454b9ad928SBarry Smith { 1464b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 147dfbe8321SBarry Smith PetscErrorCode ierr; 1489dcbbd2bSBarry Smith PetscInt nmax = 8,i; 1494b9ad928SBarry Smith PC_CompositeLink next; 150e5999256SBarry Smith char *pcs[8]; 1514b9ad928SBarry Smith PetscTruth flg; 1524b9ad928SBarry Smith 1534b9ad928SBarry Smith PetscFunctionBegin; 1544b9ad928SBarry Smith ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 1559dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_composite_type","Type of composition","PCCompositeSetType",PCCompositeTypes,(PetscEnum)jac->type,(PetscEnum*)&jac->type,&flg);CHKERRQ(ierr); 1564b9ad928SBarry Smith ierr = PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);CHKERRQ(ierr); 1574b9ad928SBarry Smith if (flg) { 1584b9ad928SBarry Smith ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr); 1594b9ad928SBarry Smith } 1604b9ad928SBarry Smith ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 1614b9ad928SBarry Smith if (flg) { 1624b9ad928SBarry Smith for (i=0; i<nmax; i++) { 1634b9ad928SBarry Smith ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 1644b9ad928SBarry Smith } 1654b9ad928SBarry Smith } 1664b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 1674b9ad928SBarry Smith 1684b9ad928SBarry Smith next = jac->head; 1694b9ad928SBarry Smith while (next) { 1704b9ad928SBarry Smith ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 1714b9ad928SBarry Smith next = next->next; 1724b9ad928SBarry Smith } 1734b9ad928SBarry Smith PetscFunctionReturn(0); 1744b9ad928SBarry Smith } 1754b9ad928SBarry Smith 1764b9ad928SBarry Smith #undef __FUNCT__ 1774b9ad928SBarry Smith #define __FUNCT__ "PCView_Composite" 1786849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer) 1794b9ad928SBarry Smith { 1804b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 181dfbe8321SBarry Smith PetscErrorCode ierr; 1824b9ad928SBarry Smith PC_CompositeLink next = jac->head; 18332077d6dSBarry Smith PetscTruth iascii; 1844b9ad928SBarry Smith 1854b9ad928SBarry Smith PetscFunctionBegin; 18632077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 18732077d6dSBarry Smith if (iascii) { 1889dcbbd2bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Composite PC type - %s\n",PCCompositeTypes[jac->type]);CHKERRQ(ierr); 1894b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 1904b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 1914b9ad928SBarry Smith } else { 19279a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name); 1934b9ad928SBarry Smith } 19432077d6dSBarry Smith if (iascii) { 1954b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 1964b9ad928SBarry Smith } 1974b9ad928SBarry Smith while (next) { 1984b9ad928SBarry Smith ierr = PCView(next->pc,viewer);CHKERRQ(ierr); 1994b9ad928SBarry Smith next = next->next; 2004b9ad928SBarry Smith } 20132077d6dSBarry Smith if (iascii) { 2024b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 2034b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 2044b9ad928SBarry Smith } 2054b9ad928SBarry Smith PetscFunctionReturn(0); 2064b9ad928SBarry Smith } 2074b9ad928SBarry Smith 2084b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 2094b9ad928SBarry Smith 2104b9ad928SBarry Smith EXTERN_C_BEGIN 2114b9ad928SBarry Smith #undef __FUNCT__ 2124b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite" 213dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 2144b9ad928SBarry Smith { 2154b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 2164b9ad928SBarry Smith PetscFunctionBegin; 2174b9ad928SBarry Smith jac->alpha = alpha; 2184b9ad928SBarry Smith PetscFunctionReturn(0); 2194b9ad928SBarry Smith } 2204b9ad928SBarry Smith EXTERN_C_END 2214b9ad928SBarry Smith 2224b9ad928SBarry Smith EXTERN_C_BEGIN 2234b9ad928SBarry Smith #undef __FUNCT__ 2244b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType_Composite" 225dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType_Composite(PC pc,PCCompositeType type) 2264b9ad928SBarry Smith { 2274b9ad928SBarry Smith PetscFunctionBegin; 2284b9ad928SBarry Smith if (type == PC_COMPOSITE_ADDITIVE) { 2294b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 2304b9ad928SBarry Smith } else if (type == PC_COMPOSITE_MULTIPLICATIVE) { 2314b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Multiplicative; 2324b9ad928SBarry Smith } else if (type == PC_COMPOSITE_SPECIAL) { 2334b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Special; 2344b9ad928SBarry Smith } else { 2351302d50aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Unkown composite preconditioner type"); 2364b9ad928SBarry Smith } 2374b9ad928SBarry Smith PetscFunctionReturn(0); 2384b9ad928SBarry Smith } 2394b9ad928SBarry Smith EXTERN_C_END 2404b9ad928SBarry Smith 2414b9ad928SBarry Smith EXTERN_C_BEGIN 2424b9ad928SBarry Smith #undef __FUNCT__ 2434b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC_Composite" 244dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC_Composite(PC pc,PCType type) 2454b9ad928SBarry Smith { 2464b9ad928SBarry Smith PC_Composite *jac; 2475a9f2f41SSatish Balay PC_CompositeLink next,ilink; 248dfbe8321SBarry Smith PetscErrorCode ierr; 24979416396SBarry Smith PetscInt cnt = 0; 2502dcb1b2aSMatthew Knepley const char *prefix; 2512dcb1b2aSMatthew Knepley char newprefix[8]; 2524b9ad928SBarry Smith 2534b9ad928SBarry Smith PetscFunctionBegin; 2545a9f2f41SSatish Balay ierr = PetscNew(struct _PC_CompositeLink,&ilink);CHKERRQ(ierr); 2555a9f2f41SSatish Balay ilink->next = 0; 2565a9f2f41SSatish Balay ierr = PCCreate(pc->comm,&ilink->pc);CHKERRQ(ierr); 2574b9ad928SBarry Smith 2584b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 2594b9ad928SBarry Smith next = jac->head; 2604b9ad928SBarry Smith if (!next) { 2615a9f2f41SSatish Balay jac->head = ilink; 2624b9ad928SBarry Smith } else { 2634b9ad928SBarry Smith cnt++; 2644b9ad928SBarry Smith while (next->next) { 2654b9ad928SBarry Smith next = next->next; 2664b9ad928SBarry Smith cnt++; 2674b9ad928SBarry Smith } 2685a9f2f41SSatish Balay next->next = ilink; 2694b9ad928SBarry Smith } 2704b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 2715a9f2f41SSatish Balay ierr = PCSetOptionsPrefix(ilink->pc,prefix);CHKERRQ(ierr); 27213f74950SBarry Smith sprintf(newprefix,"sub_%d_",(int)cnt); 2735a9f2f41SSatish Balay ierr = PCAppendOptionsPrefix(ilink->pc,newprefix);CHKERRQ(ierr); 2744b9ad928SBarry Smith /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 2755a9f2f41SSatish Balay ierr = PCSetType(ilink->pc,type);CHKERRQ(ierr); 2764b9ad928SBarry Smith 2774b9ad928SBarry Smith PetscFunctionReturn(0); 2784b9ad928SBarry Smith } 2794b9ad928SBarry Smith EXTERN_C_END 2804b9ad928SBarry Smith 2814b9ad928SBarry Smith EXTERN_C_BEGIN 2824b9ad928SBarry Smith #undef __FUNCT__ 2834b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC_Composite" 284dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC_Composite(PC pc,PetscInt n,PC *subpc) 2854b9ad928SBarry Smith { 2864b9ad928SBarry Smith PC_Composite *jac; 2874b9ad928SBarry Smith PC_CompositeLink next; 28879416396SBarry Smith PetscInt i; 2894b9ad928SBarry Smith 2904b9ad928SBarry Smith PetscFunctionBegin; 2914b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 2924b9ad928SBarry Smith next = jac->head; 2934b9ad928SBarry Smith for (i=0; i<n; i++) { 2944b9ad928SBarry Smith if (!next->next) { 2951302d50aSBarry Smith SETERRQ(PETSC_ERR_ARG_INCOMP,"Not enough PCs in composite preconditioner"); 2964b9ad928SBarry Smith } 2974b9ad928SBarry Smith next = next->next; 2984b9ad928SBarry Smith } 2994b9ad928SBarry Smith *subpc = next->pc; 3004b9ad928SBarry Smith PetscFunctionReturn(0); 3014b9ad928SBarry Smith } 3024b9ad928SBarry Smith EXTERN_C_END 3034b9ad928SBarry Smith 3044b9ad928SBarry Smith EXTERN_C_BEGIN 3054b9ad928SBarry Smith #undef __FUNCT__ 3064b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue_Composite" 307dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue_Composite(PC pc) 3084b9ad928SBarry Smith { 3094b9ad928SBarry Smith PC_Composite *jac; 3104b9ad928SBarry Smith 3114b9ad928SBarry Smith PetscFunctionBegin; 3124b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 3134b9ad928SBarry Smith jac->use_true_matrix = PETSC_TRUE; 3144b9ad928SBarry Smith PetscFunctionReturn(0); 3154b9ad928SBarry Smith } 3164b9ad928SBarry Smith EXTERN_C_END 3174b9ad928SBarry Smith 3184b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */ 3194b9ad928SBarry Smith #undef __FUNCT__ 3204b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType" 321*f39d8e23SSatish Balay /*@ 3224b9ad928SBarry Smith PCCompositeSetType - Sets the type of composite preconditioner. 3234b9ad928SBarry Smith 3244b9ad928SBarry Smith Collective on PC 3254b9ad928SBarry Smith 3264b9ad928SBarry Smith Input Parameter: 3274b9ad928SBarry Smith . pc - the preconditioner context 3284b9ad928SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 3294b9ad928SBarry Smith 3304b9ad928SBarry Smith Options Database Key: 3314b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 3324b9ad928SBarry Smith 3334b9ad928SBarry Smith Level: Developer 3344b9ad928SBarry Smith 3354b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 3364b9ad928SBarry Smith @*/ 337dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetType(PC pc,PCCompositeType type) 3384b9ad928SBarry Smith { 339dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,PCCompositeType); 3404b9ad928SBarry Smith 3414b9ad928SBarry Smith PetscFunctionBegin; 3424482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3434b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 3444b9ad928SBarry Smith if (f) { 3454b9ad928SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 3464b9ad928SBarry Smith } 3474b9ad928SBarry Smith PetscFunctionReturn(0); 3484b9ad928SBarry Smith } 3494b9ad928SBarry Smith 3504b9ad928SBarry Smith #undef __FUNCT__ 3514b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha" 352*f39d8e23SSatish Balay /*@ 3534b9ad928SBarry Smith PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner 3544b9ad928SBarry Smith for alphaI + R + S 3554b9ad928SBarry Smith 3564b9ad928SBarry Smith Collective on PC 3574b9ad928SBarry Smith 3584b9ad928SBarry Smith Input Parameter: 3594b9ad928SBarry Smith + pc - the preconditioner context 3604b9ad928SBarry Smith - alpha - scale on identity 3614b9ad928SBarry Smith 3624b9ad928SBarry Smith Level: Developer 3634b9ad928SBarry Smith 3644b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 3654b9ad928SBarry Smith @*/ 366dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 3674b9ad928SBarry Smith { 368dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscScalar); 3694b9ad928SBarry Smith 3704b9ad928SBarry Smith PetscFunctionBegin; 3714482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3724b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);CHKERRQ(ierr); 3734b9ad928SBarry Smith if (f) { 3744b9ad928SBarry Smith ierr = (*f)(pc,alpha);CHKERRQ(ierr); 3754b9ad928SBarry Smith } 3764b9ad928SBarry Smith PetscFunctionReturn(0); 3774b9ad928SBarry Smith } 3784b9ad928SBarry Smith 3794b9ad928SBarry Smith #undef __FUNCT__ 3804b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC" 3814b9ad928SBarry Smith /*@C 3824b9ad928SBarry Smith PCCompositeAddPC - Adds another PC to the composite PC. 3834b9ad928SBarry Smith 3844b9ad928SBarry Smith Collective on PC 3854b9ad928SBarry Smith 3864b9ad928SBarry Smith Input Parameters: 3874b9ad928SBarry Smith . pc - the preconditioner context 3884b9ad928SBarry Smith . type - the type of the new preconditioner 3894b9ad928SBarry Smith 3904b9ad928SBarry Smith Level: Developer 3914b9ad928SBarry Smith 3924b9ad928SBarry Smith .keywords: PC, composite preconditioner, add 3934b9ad928SBarry Smith @*/ 394dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeAddPC(PC pc,PCType type) 3954b9ad928SBarry Smith { 396dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,PCType); 3974b9ad928SBarry Smith 3984b9ad928SBarry Smith PetscFunctionBegin; 3994482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4004b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);CHKERRQ(ierr); 4014b9ad928SBarry Smith if (f) { 4024b9ad928SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 4034b9ad928SBarry Smith } 4044b9ad928SBarry Smith PetscFunctionReturn(0); 4054b9ad928SBarry Smith } 4064b9ad928SBarry Smith 4074b9ad928SBarry Smith #undef __FUNCT__ 4084b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC" 409*f39d8e23SSatish Balay /*@ 4104b9ad928SBarry Smith PCCompositeGetPC - Gets one of the PC objects in the composite PC. 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith Not Collective 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith Input Parameter: 4154b9ad928SBarry Smith . pc - the preconditioner context 4164b9ad928SBarry Smith . n - the number of the pc requested 4174b9ad928SBarry Smith 4184b9ad928SBarry Smith Output Parameters: 4194b9ad928SBarry Smith . subpc - the PC requested 4204b9ad928SBarry Smith 4214b9ad928SBarry Smith Level: Developer 4224b9ad928SBarry Smith 4234b9ad928SBarry Smith .keywords: PC, get, composite preconditioner, sub preconditioner 4244b9ad928SBarry Smith 4254b9ad928SBarry Smith .seealso: PCCompositeAddPC() 4264b9ad928SBarry Smith @*/ 427dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeGetPC(PC pc,PetscInt n,PC *subpc) 4284b9ad928SBarry Smith { 42913f74950SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PC *); 4304b9ad928SBarry Smith 4314b9ad928SBarry Smith PetscFunctionBegin; 4324482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4334482741eSBarry Smith PetscValidPointer(subpc,3); 4344b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 4354b9ad928SBarry Smith if (f) { 4364b9ad928SBarry Smith ierr = (*f)(pc,n,subpc);CHKERRQ(ierr); 4374b9ad928SBarry Smith } else { 4381302d50aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get pc, not composite type"); 4394b9ad928SBarry Smith } 4404b9ad928SBarry Smith PetscFunctionReturn(0); 4414b9ad928SBarry Smith } 4424b9ad928SBarry Smith 4434b9ad928SBarry Smith #undef __FUNCT__ 4444b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue" 4454b9ad928SBarry Smith /*@ 4464b9ad928SBarry Smith PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than 4474b9ad928SBarry Smith the matrix used to define the preconditioner) is used to compute 4484b9ad928SBarry Smith the residual when the multiplicative scheme is used. 4494b9ad928SBarry Smith 4504b9ad928SBarry Smith Collective on PC 4514b9ad928SBarry Smith 4524b9ad928SBarry Smith Input Parameters: 4534b9ad928SBarry Smith . pc - the preconditioner context 4544b9ad928SBarry Smith 4554b9ad928SBarry Smith Options Database Key: 4564b9ad928SBarry Smith . -pc_composite_true - Activates PCCompositeSetUseTrue() 4574b9ad928SBarry Smith 4584b9ad928SBarry Smith Note: 4594b9ad928SBarry Smith For the common case in which the preconditioning and linear 4604b9ad928SBarry Smith system matrices are identical, this routine is unnecessary. 4614b9ad928SBarry Smith 4624b9ad928SBarry Smith Level: Developer 4634b9ad928SBarry Smith 4644b9ad928SBarry Smith .keywords: PC, composite preconditioner, set, true, flag 4654b9ad928SBarry Smith 4664b9ad928SBarry Smith .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue() 4674b9ad928SBarry Smith @*/ 468dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCompositeSetUseTrue(PC pc) 4694b9ad928SBarry Smith { 470dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC); 4714b9ad928SBarry Smith 4724b9ad928SBarry Smith PetscFunctionBegin; 4734482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4744b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr); 4754b9ad928SBarry Smith if (f) { 4764b9ad928SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 4774b9ad928SBarry Smith } 4784b9ad928SBarry Smith PetscFunctionReturn(0); 4794b9ad928SBarry Smith } 4804b9ad928SBarry Smith 4814b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/ 4824b9ad928SBarry Smith 4834b9ad928SBarry Smith /*MC 4844b9ad928SBarry Smith PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 4854b9ad928SBarry Smith 4864b9ad928SBarry Smith Options Database Keys: 4874b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 4884b9ad928SBarry Smith . -pc_composite_true - Activates PCCompositeSetUseTrue() 4894b9ad928SBarry Smith 4904b9ad928SBarry Smith Level: intermediate 4914b9ad928SBarry Smith 4924b9ad928SBarry Smith Concepts: composing solvers 4934b9ad928SBarry Smith 4944b9ad928SBarry Smith Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 4954b9ad928SBarry Smith inner PCs to be PCKSP. 4964b9ad928SBarry Smith Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 49779416396SBarry Smith the incorrect answer) unless you use KSPFGMRES as the outter Krylov method 4984b9ad928SBarry Smith 4994b9ad928SBarry Smith 5004b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 5014b9ad928SBarry Smith PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 5024b9ad928SBarry Smith PCCompositeGetPC(), PCCompositeSetUseTrue() 5034b9ad928SBarry Smith 5044b9ad928SBarry Smith M*/ 5054b9ad928SBarry Smith 5064b9ad928SBarry Smith EXTERN_C_BEGIN 5074b9ad928SBarry Smith #undef __FUNCT__ 5084b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Composite" 509dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Composite(PC pc) 5104b9ad928SBarry Smith { 511dfbe8321SBarry Smith PetscErrorCode ierr; 5124b9ad928SBarry Smith PC_Composite *jac; 5134b9ad928SBarry Smith 5144b9ad928SBarry Smith PetscFunctionBegin; 5154b9ad928SBarry Smith ierr = PetscNew(PC_Composite,&jac);CHKERRQ(ierr); 51652e6d16bSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PC_Composite));CHKERRQ(ierr); 5174b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 5184b9ad928SBarry Smith pc->ops->setup = PCSetUp_Composite; 5194b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Composite; 5204b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Composite; 5214b9ad928SBarry Smith pc->ops->view = PCView_Composite; 5224b9ad928SBarry Smith pc->ops->applyrichardson = 0; 5234b9ad928SBarry Smith 5244b9ad928SBarry Smith pc->data = (void*)jac; 5254b9ad928SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 5264b9ad928SBarry Smith jac->work1 = 0; 5274b9ad928SBarry Smith jac->work2 = 0; 5284b9ad928SBarry Smith jac->head = 0; 5294b9ad928SBarry Smith 5304b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite", 5314b9ad928SBarry Smith PCCompositeSetType_Composite);CHKERRQ(ierr); 5324b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite", 5334b9ad928SBarry Smith PCCompositeAddPC_Composite);CHKERRQ(ierr); 5344b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite", 5354b9ad928SBarry Smith PCCompositeGetPC_Composite);CHKERRQ(ierr); 5364b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite", 5374b9ad928SBarry Smith PCCompositeSetUseTrue_Composite);CHKERRQ(ierr); 5384b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite", 5394b9ad928SBarry Smith PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr); 5404b9ad928SBarry Smith 5414b9ad928SBarry Smith PetscFunctionReturn(0); 5424b9ad928SBarry Smith } 5434b9ad928SBarry Smith EXTERN_C_END 5444b9ad928SBarry Smith 545