14b9ad928SBarry Smith /* 24b9ad928SBarry Smith Defines a preconditioner that can consist of a collection of PCs 34b9ad928SBarry Smith */ 44b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 54b9ad928SBarry Smith #include "petscksp.h" /*I "petscksp.h" I*/ 64b9ad928SBarry Smith 74b9ad928SBarry Smith typedef struct _PC_CompositeLink *PC_CompositeLink; 84b9ad928SBarry Smith struct _PC_CompositeLink { 94b9ad928SBarry Smith PC pc; 104b9ad928SBarry Smith PC_CompositeLink next; 114b9ad928SBarry Smith }; 124b9ad928SBarry Smith 134b9ad928SBarry Smith typedef struct { 144b9ad928SBarry Smith PC_CompositeLink head; 154b9ad928SBarry Smith PCCompositeType type; 164b9ad928SBarry Smith Vec work1; 174b9ad928SBarry Smith Vec work2; 184b9ad928SBarry Smith PetscScalar alpha; 194b9ad928SBarry Smith PetscTruth use_true_matrix; 204b9ad928SBarry Smith } PC_Composite; 214b9ad928SBarry Smith 224b9ad928SBarry Smith #undef __FUNCT__ 234b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Multiplicative" 24*6849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y) 254b9ad928SBarry Smith { 26dfbe8321SBarry Smith PetscErrorCode ierr; 274b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 284b9ad928SBarry Smith PC_CompositeLink next = jac->head; 294b9ad928SBarry Smith PetscScalar one = 1.0,mone = -1.0; 304b9ad928SBarry Smith Mat mat = pc->pmat; 314b9ad928SBarry Smith 324b9ad928SBarry Smith PetscFunctionBegin; 334b9ad928SBarry Smith if (!next) { 344b9ad928SBarry Smith SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 354b9ad928SBarry Smith } 364b9ad928SBarry Smith if (next->next && !jac->work2) { /* allocate second work vector */ 374b9ad928SBarry Smith ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr); 384b9ad928SBarry Smith } 39d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 404b9ad928SBarry Smith if (jac->use_true_matrix) mat = pc->mat; 414b9ad928SBarry Smith while (next->next) { 424b9ad928SBarry Smith next = next->next; 434b9ad928SBarry Smith ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr); 444b9ad928SBarry Smith ierr = VecWAXPY(&mone,jac->work1,x,jac->work2);CHKERRQ(ierr); 45d8fd42c4SBarry Smith ierr = PCApply(next->pc,jac->work2,jac->work1);CHKERRQ(ierr); 464b9ad928SBarry Smith ierr = VecAXPY(&one,jac->work1,y);CHKERRQ(ierr); 474b9ad928SBarry Smith } 484b9ad928SBarry Smith 494b9ad928SBarry Smith PetscFunctionReturn(0); 504b9ad928SBarry Smith } 514b9ad928SBarry Smith 524b9ad928SBarry Smith /* 534b9ad928SBarry Smith This is very special for a matrix of the form alpha I + R + S 544b9ad928SBarry Smith where first preconditioner is built from alpha I + S and second from 554b9ad928SBarry Smith alpha I + R 564b9ad928SBarry Smith */ 574b9ad928SBarry Smith #undef __FUNCT__ 584b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Special" 59*6849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Special(PC pc,Vec x,Vec y) 604b9ad928SBarry Smith { 61dfbe8321SBarry Smith PetscErrorCode ierr; 624b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 634b9ad928SBarry Smith PC_CompositeLink next = jac->head; 644b9ad928SBarry Smith 654b9ad928SBarry Smith PetscFunctionBegin; 664b9ad928SBarry Smith if (!next) { 674b9ad928SBarry Smith SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 684b9ad928SBarry Smith } 694b9ad928SBarry Smith if (!next->next || next->next->next) { 704b9ad928SBarry Smith SETERRQ(1,"Special composite preconditioners requires exactly two PCs"); 714b9ad928SBarry Smith } 724b9ad928SBarry Smith 73d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 74d8fd42c4SBarry Smith ierr = PCApply(next->next->pc,jac->work1,y);CHKERRQ(ierr); 754b9ad928SBarry Smith PetscFunctionReturn(0); 764b9ad928SBarry Smith } 774b9ad928SBarry Smith 784b9ad928SBarry Smith #undef __FUNCT__ 794b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Additive" 80*6849ba73SBarry Smith static PetscErrorCode PCApply_Composite_Additive(PC pc,Vec x,Vec y) 814b9ad928SBarry Smith { 82dfbe8321SBarry Smith PetscErrorCode ierr; 834b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 844b9ad928SBarry Smith PC_CompositeLink next = jac->head; 854b9ad928SBarry Smith PetscScalar one = 1.0; 864b9ad928SBarry Smith 874b9ad928SBarry Smith PetscFunctionBegin; 884b9ad928SBarry Smith if (!next) { 894b9ad928SBarry Smith SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 904b9ad928SBarry Smith } 91d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,y);CHKERRQ(ierr); 924b9ad928SBarry Smith while (next->next) { 934b9ad928SBarry Smith next = next->next; 94d8fd42c4SBarry Smith ierr = PCApply(next->pc,x,jac->work1);CHKERRQ(ierr); 954b9ad928SBarry Smith ierr = VecAXPY(&one,jac->work1,y);CHKERRQ(ierr); 964b9ad928SBarry Smith } 974b9ad928SBarry Smith PetscFunctionReturn(0); 984b9ad928SBarry Smith } 994b9ad928SBarry Smith 1004b9ad928SBarry Smith #undef __FUNCT__ 1014b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Composite" 102*6849ba73SBarry Smith static PetscErrorCode PCSetUp_Composite(PC pc) 1034b9ad928SBarry Smith { 104dfbe8321SBarry Smith PetscErrorCode ierr; 1054b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 1064b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1074b9ad928SBarry Smith 1084b9ad928SBarry Smith PetscFunctionBegin; 1094b9ad928SBarry Smith if (!jac->work1) { 11023ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&jac->work1,0);CHKERRQ(ierr); 1114b9ad928SBarry Smith } 1124b9ad928SBarry Smith while (next) { 1134b9ad928SBarry Smith ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1144b9ad928SBarry Smith next = next->next; 1154b9ad928SBarry Smith } 1164b9ad928SBarry Smith 1174b9ad928SBarry Smith PetscFunctionReturn(0); 1184b9ad928SBarry Smith } 1194b9ad928SBarry Smith 1204b9ad928SBarry Smith #undef __FUNCT__ 1214b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Composite" 122*6849ba73SBarry Smith static PetscErrorCode PCDestroy_Composite(PC pc) 1234b9ad928SBarry Smith { 1244b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 125dfbe8321SBarry Smith PetscErrorCode ierr; 1264b9ad928SBarry Smith PC_CompositeLink next = jac->head; 1274b9ad928SBarry Smith 1284b9ad928SBarry Smith PetscFunctionBegin; 1294b9ad928SBarry Smith while (next) { 1304b9ad928SBarry Smith ierr = PCDestroy(next->pc);CHKERRQ(ierr); 1314b9ad928SBarry Smith next = next->next; 1324b9ad928SBarry Smith } 1334b9ad928SBarry Smith 1344b9ad928SBarry Smith if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);} 1354b9ad928SBarry Smith if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);} 1364b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 1374b9ad928SBarry Smith PetscFunctionReturn(0); 1384b9ad928SBarry Smith } 1394b9ad928SBarry Smith 1404b9ad928SBarry Smith #undef __FUNCT__ 1414b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite" 142*6849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Composite(PC pc) 1434b9ad928SBarry Smith { 1444b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 145dfbe8321SBarry Smith PetscErrorCode ierr; 146dfbe8321SBarry Smith int 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" 180*6849ba73SBarry Smith static PetscErrorCode PCView_Composite(PC pc,PetscViewer viewer) 1814b9ad928SBarry Smith { 1824b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 183dfbe8321SBarry Smith PetscErrorCode ierr; 1844b9ad928SBarry Smith PC_CompositeLink next = jac->head; 18532077d6dSBarry Smith PetscTruth iascii; 1864b9ad928SBarry Smith 1874b9ad928SBarry Smith PetscFunctionBegin; 18832077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 18932077d6dSBarry Smith if (iascii) { 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 } 19532077d6dSBarry Smith if (iascii) { 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 } 20232077d6dSBarry Smith if (iascii) { 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" 214dfbe8321SBarry Smith PetscErrorCode 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" 226dfbe8321SBarry Smith PetscErrorCode 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" 245dfbe8321SBarry Smith PetscErrorCode PCCompositeAddPC_Composite(PC pc,PCType type) 2464b9ad928SBarry Smith { 2474b9ad928SBarry Smith PC_Composite *jac; 2484b9ad928SBarry Smith PC_CompositeLink next,link; 249dfbe8321SBarry Smith PetscErrorCode ierr; 250dfbe8321SBarry Smith int cnt = 0; 2514b9ad928SBarry Smith char *prefix,newprefix[8]; 2524b9ad928SBarry Smith 2534b9ad928SBarry Smith PetscFunctionBegin; 2544b9ad928SBarry Smith ierr = PetscNew(struct _PC_CompositeLink,&link);CHKERRQ(ierr); 2554b9ad928SBarry Smith link->next = 0; 2564b9ad928SBarry Smith ierr = PCCreate(pc->comm,&link->pc);CHKERRQ(ierr); 2574b9ad928SBarry Smith 2584b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 2594b9ad928SBarry Smith next = jac->head; 2604b9ad928SBarry Smith if (!next) { 2614b9ad928SBarry Smith jac->head = link; 2624b9ad928SBarry Smith } else { 2634b9ad928SBarry Smith cnt++; 2644b9ad928SBarry Smith while (next->next) { 2654b9ad928SBarry Smith next = next->next; 2664b9ad928SBarry Smith cnt++; 2674b9ad928SBarry Smith } 2684b9ad928SBarry Smith next->next = link; 2694b9ad928SBarry Smith } 2704b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 2714b9ad928SBarry Smith ierr = PCSetOptionsPrefix(link->pc,prefix);CHKERRQ(ierr); 2724b9ad928SBarry Smith sprintf(newprefix,"sub_%d_",cnt); 2734b9ad928SBarry Smith ierr = PCAppendOptionsPrefix(link->pc,newprefix);CHKERRQ(ierr); 2744b9ad928SBarry Smith /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 2754b9ad928SBarry Smith ierr = PCSetType(link->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" 284dfbe8321SBarry Smith PetscErrorCode PCCompositeGetPC_Composite(PC pc,int n,PC *subpc) 2854b9ad928SBarry Smith { 2864b9ad928SBarry Smith PC_Composite *jac; 2874b9ad928SBarry Smith PC_CompositeLink next; 2884b9ad928SBarry Smith int 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) { 2954b9ad928SBarry Smith SETERRQ(1,"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" 307dfbe8321SBarry Smith PetscErrorCode 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" 3214b9ad928SBarry Smith /*@C 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 @*/ 337dfbe8321SBarry Smith PetscErrorCode 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" 3524b9ad928SBarry Smith /*@C 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 @*/ 366dfbe8321SBarry Smith PetscErrorCode 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 @*/ 394dfbe8321SBarry Smith PetscErrorCode 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" 4094b9ad928SBarry Smith /*@C 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 @*/ 427dfbe8321SBarry Smith PetscErrorCode PCCompositeGetPC(PC pc,int n,PC *subpc) 4284b9ad928SBarry Smith { 429dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,int,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 { 4384b9ad928SBarry Smith SETERRQ(1,"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 @*/ 468dfbe8321SBarry Smith PetscErrorCode 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 4974b9ad928SBarry Smith the incorrect answer) unless you use KSPFGMRES as the other 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" 509dfbe8321SBarry Smith PetscErrorCode 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); 5164b9ad928SBarry Smith PetscLogObjectMemory(pc,sizeof(PC_Composite)); 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