1*4b9ad928SBarry Smith /*$Id: composite.c,v 1.45 2001/08/07 03:03:39 balay Exp $*/ 2*4b9ad928SBarry Smith /* 3*4b9ad928SBarry Smith Defines a preconditioner that can consist of a collection of PCs 4*4b9ad928SBarry Smith */ 5*4b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 6*4b9ad928SBarry Smith #include "petscksp.h" /*I "petscksp.h" I*/ 7*4b9ad928SBarry Smith 8*4b9ad928SBarry Smith typedef struct _PC_CompositeLink *PC_CompositeLink; 9*4b9ad928SBarry Smith struct _PC_CompositeLink { 10*4b9ad928SBarry Smith PC pc; 11*4b9ad928SBarry Smith PC_CompositeLink next; 12*4b9ad928SBarry Smith }; 13*4b9ad928SBarry Smith 14*4b9ad928SBarry Smith typedef struct { 15*4b9ad928SBarry Smith PC_CompositeLink head; 16*4b9ad928SBarry Smith PCCompositeType type; 17*4b9ad928SBarry Smith Vec work1; 18*4b9ad928SBarry Smith Vec work2; 19*4b9ad928SBarry Smith PetscScalar alpha; 20*4b9ad928SBarry Smith PetscTruth use_true_matrix; 21*4b9ad928SBarry Smith } PC_Composite; 22*4b9ad928SBarry Smith 23*4b9ad928SBarry Smith #undef __FUNCT__ 24*4b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Multiplicative" 25*4b9ad928SBarry Smith static int PCApply_Composite_Multiplicative(PC pc,Vec x,Vec y) 26*4b9ad928SBarry Smith { 27*4b9ad928SBarry Smith int ierr; 28*4b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 29*4b9ad928SBarry Smith PC_CompositeLink next = jac->head; 30*4b9ad928SBarry Smith PetscScalar one = 1.0,mone = -1.0; 31*4b9ad928SBarry Smith Mat mat = pc->pmat; 32*4b9ad928SBarry Smith 33*4b9ad928SBarry Smith PetscFunctionBegin; 34*4b9ad928SBarry Smith if (!next) { 35*4b9ad928SBarry Smith SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 36*4b9ad928SBarry Smith } 37*4b9ad928SBarry Smith if (next->next && !jac->work2) { /* allocate second work vector */ 38*4b9ad928SBarry Smith ierr = VecDuplicate(jac->work1,&jac->work2);CHKERRQ(ierr); 39*4b9ad928SBarry Smith } 40*4b9ad928SBarry Smith ierr = PCApply(next->pc,x,y,PC_LEFT);CHKERRQ(ierr); 41*4b9ad928SBarry Smith if (jac->use_true_matrix) mat = pc->mat; 42*4b9ad928SBarry Smith while (next->next) { 43*4b9ad928SBarry Smith next = next->next; 44*4b9ad928SBarry Smith ierr = MatMult(mat,y,jac->work1);CHKERRQ(ierr); 45*4b9ad928SBarry Smith ierr = VecWAXPY(&mone,jac->work1,x,jac->work2);CHKERRQ(ierr); 46*4b9ad928SBarry Smith ierr = PCApply(next->pc,jac->work2,jac->work1,PC_LEFT);CHKERRQ(ierr); 47*4b9ad928SBarry Smith ierr = VecAXPY(&one,jac->work1,y);CHKERRQ(ierr); 48*4b9ad928SBarry Smith } 49*4b9ad928SBarry Smith 50*4b9ad928SBarry Smith PetscFunctionReturn(0); 51*4b9ad928SBarry Smith } 52*4b9ad928SBarry Smith 53*4b9ad928SBarry Smith /* 54*4b9ad928SBarry Smith This is very special for a matrix of the form alpha I + R + S 55*4b9ad928SBarry Smith where first preconditioner is built from alpha I + S and second from 56*4b9ad928SBarry Smith alpha I + R 57*4b9ad928SBarry Smith */ 58*4b9ad928SBarry Smith #undef __FUNCT__ 59*4b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Special" 60*4b9ad928SBarry Smith static int PCApply_Composite_Special(PC pc,Vec x,Vec y) 61*4b9ad928SBarry Smith { 62*4b9ad928SBarry Smith int ierr; 63*4b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 64*4b9ad928SBarry Smith PC_CompositeLink next = jac->head; 65*4b9ad928SBarry Smith 66*4b9ad928SBarry Smith PetscFunctionBegin; 67*4b9ad928SBarry Smith if (!next) { 68*4b9ad928SBarry Smith SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 69*4b9ad928SBarry Smith } 70*4b9ad928SBarry Smith if (!next->next || next->next->next) { 71*4b9ad928SBarry Smith SETERRQ(1,"Special composite preconditioners requires exactly two PCs"); 72*4b9ad928SBarry Smith } 73*4b9ad928SBarry Smith 74*4b9ad928SBarry Smith ierr = PCApply(next->pc,x,jac->work1,PC_LEFT);CHKERRQ(ierr); 75*4b9ad928SBarry Smith ierr = PCApply(next->next->pc,jac->work1,y,PC_LEFT);CHKERRQ(ierr); 76*4b9ad928SBarry Smith PetscFunctionReturn(0); 77*4b9ad928SBarry Smith } 78*4b9ad928SBarry Smith 79*4b9ad928SBarry Smith #undef __FUNCT__ 80*4b9ad928SBarry Smith #define __FUNCT__ "PCApply_Composite_Additive" 81*4b9ad928SBarry Smith static int PCApply_Composite_Additive(PC pc,Vec x,Vec y) 82*4b9ad928SBarry Smith { 83*4b9ad928SBarry Smith int ierr; 84*4b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 85*4b9ad928SBarry Smith PC_CompositeLink next = jac->head; 86*4b9ad928SBarry Smith PetscScalar one = 1.0; 87*4b9ad928SBarry Smith 88*4b9ad928SBarry Smith PetscFunctionBegin; 89*4b9ad928SBarry Smith if (!next) { 90*4b9ad928SBarry Smith SETERRQ(1,"No composite preconditioners supplied via PCCompositeAddPC()"); 91*4b9ad928SBarry Smith } 92*4b9ad928SBarry Smith ierr = PCApply(next->pc,x,y,PC_LEFT);CHKERRQ(ierr); 93*4b9ad928SBarry Smith while (next->next) { 94*4b9ad928SBarry Smith next = next->next; 95*4b9ad928SBarry Smith ierr = PCApply(next->pc,x,jac->work1,PC_LEFT);CHKERRQ(ierr); 96*4b9ad928SBarry Smith ierr = VecAXPY(&one,jac->work1,y);CHKERRQ(ierr); 97*4b9ad928SBarry Smith } 98*4b9ad928SBarry Smith PetscFunctionReturn(0); 99*4b9ad928SBarry Smith } 100*4b9ad928SBarry Smith 101*4b9ad928SBarry Smith #undef __FUNCT__ 102*4b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Composite" 103*4b9ad928SBarry Smith static int PCSetUp_Composite(PC pc) 104*4b9ad928SBarry Smith { 105*4b9ad928SBarry Smith int ierr; 106*4b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 107*4b9ad928SBarry Smith PC_CompositeLink next = jac->head; 108*4b9ad928SBarry Smith 109*4b9ad928SBarry Smith PetscFunctionBegin; 110*4b9ad928SBarry Smith if (!jac->work1) { 111*4b9ad928SBarry Smith ierr = VecDuplicate(pc->vec,&jac->work1);CHKERRQ(ierr); 112*4b9ad928SBarry Smith } 113*4b9ad928SBarry Smith while (next) { 114*4b9ad928SBarry Smith ierr = PCSetOperators(next->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 115*4b9ad928SBarry Smith ierr = PCSetVector(next->pc,jac->work1);CHKERRQ(ierr); 116*4b9ad928SBarry Smith next = next->next; 117*4b9ad928SBarry Smith } 118*4b9ad928SBarry Smith 119*4b9ad928SBarry Smith PetscFunctionReturn(0); 120*4b9ad928SBarry Smith } 121*4b9ad928SBarry Smith 122*4b9ad928SBarry Smith #undef __FUNCT__ 123*4b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Composite" 124*4b9ad928SBarry Smith static int PCDestroy_Composite(PC pc) 125*4b9ad928SBarry Smith { 126*4b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 127*4b9ad928SBarry Smith int ierr; 128*4b9ad928SBarry Smith PC_CompositeLink next = jac->head; 129*4b9ad928SBarry Smith 130*4b9ad928SBarry Smith PetscFunctionBegin; 131*4b9ad928SBarry Smith while (next) { 132*4b9ad928SBarry Smith ierr = PCDestroy(next->pc);CHKERRQ(ierr); 133*4b9ad928SBarry Smith next = next->next; 134*4b9ad928SBarry Smith } 135*4b9ad928SBarry Smith 136*4b9ad928SBarry Smith if (jac->work1) {ierr = VecDestroy(jac->work1);CHKERRQ(ierr);} 137*4b9ad928SBarry Smith if (jac->work2) {ierr = VecDestroy(jac->work2);CHKERRQ(ierr);} 138*4b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 139*4b9ad928SBarry Smith PetscFunctionReturn(0); 140*4b9ad928SBarry Smith } 141*4b9ad928SBarry Smith 142*4b9ad928SBarry Smith #undef __FUNCT__ 143*4b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Composite" 144*4b9ad928SBarry Smith static int PCSetFromOptions_Composite(PC pc) 145*4b9ad928SBarry Smith { 146*4b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 147*4b9ad928SBarry Smith int ierr,nmax = 8,i,indx; 148*4b9ad928SBarry Smith PC_CompositeLink next; 149*4b9ad928SBarry Smith char *pcs[8],*types[] = {"multiplicative","additive","special"}; 150*4b9ad928SBarry Smith PetscTruth flg; 151*4b9ad928SBarry Smith 152*4b9ad928SBarry Smith PetscFunctionBegin; 153*4b9ad928SBarry Smith ierr = PetscOptionsHead("Composite preconditioner options");CHKERRQ(ierr); 154*4b9ad928SBarry Smith ierr = PetscOptionsEList("-pc_composite_type","Type of composition","PCCompositeSetType",types,3,types[0],&indx,&flg);CHKERRQ(ierr); 155*4b9ad928SBarry Smith if (flg) { 156*4b9ad928SBarry Smith ierr = PCCompositeSetType(pc,(PCCompositeType)indx);CHKERRQ(ierr); 157*4b9ad928SBarry Smith } 158*4b9ad928SBarry Smith ierr = PetscOptionsName("-pc_composite_true","Use true matrix for inner solves","PCCompositeSetUseTrue",&flg);CHKERRQ(ierr); 159*4b9ad928SBarry Smith if (flg) { 160*4b9ad928SBarry Smith ierr = PCCompositeSetUseTrue(pc);CHKERRQ(ierr); 161*4b9ad928SBarry Smith } 162*4b9ad928SBarry Smith ierr = PetscOptionsStringArray("-pc_composite_pcs","List of composite solvers","PCCompositeAddPC",pcs,&nmax,&flg);CHKERRQ(ierr); 163*4b9ad928SBarry Smith if (flg) { 164*4b9ad928SBarry Smith for (i=0; i<nmax; i++) { 165*4b9ad928SBarry Smith ierr = PCCompositeAddPC(pc,pcs[i]);CHKERRQ(ierr); 166*4b9ad928SBarry Smith } 167*4b9ad928SBarry Smith } 168*4b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 169*4b9ad928SBarry Smith 170*4b9ad928SBarry Smith next = jac->head; 171*4b9ad928SBarry Smith while (next) { 172*4b9ad928SBarry Smith ierr = PCSetFromOptions(next->pc);CHKERRQ(ierr); 173*4b9ad928SBarry Smith next = next->next; 174*4b9ad928SBarry Smith } 175*4b9ad928SBarry Smith PetscFunctionReturn(0); 176*4b9ad928SBarry Smith } 177*4b9ad928SBarry Smith 178*4b9ad928SBarry Smith #undef __FUNCT__ 179*4b9ad928SBarry Smith #define __FUNCT__ "PCView_Composite" 180*4b9ad928SBarry Smith static int PCView_Composite(PC pc,PetscViewer viewer) 181*4b9ad928SBarry Smith { 182*4b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 183*4b9ad928SBarry Smith int ierr; 184*4b9ad928SBarry Smith PC_CompositeLink next = jac->head; 185*4b9ad928SBarry Smith PetscTruth isascii; 186*4b9ad928SBarry Smith 187*4b9ad928SBarry Smith PetscFunctionBegin; 188*4b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 189*4b9ad928SBarry Smith if (isascii) { 190*4b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"PCs on composite preconditioner follow\n");CHKERRQ(ierr); 191*4b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 192*4b9ad928SBarry Smith } else { 193*4b9ad928SBarry Smith SETERRQ1(1,"Viewer type %s not supported for PCComposite",((PetscObject)viewer)->type_name); 194*4b9ad928SBarry Smith } 195*4b9ad928SBarry Smith if (isascii) { 196*4b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 197*4b9ad928SBarry Smith } 198*4b9ad928SBarry Smith while (next) { 199*4b9ad928SBarry Smith ierr = PCView(next->pc,viewer);CHKERRQ(ierr); 200*4b9ad928SBarry Smith next = next->next; 201*4b9ad928SBarry Smith } 202*4b9ad928SBarry Smith if (isascii) { 203*4b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 204*4b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 205*4b9ad928SBarry Smith } 206*4b9ad928SBarry Smith PetscFunctionReturn(0); 207*4b9ad928SBarry Smith } 208*4b9ad928SBarry Smith 209*4b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 210*4b9ad928SBarry Smith 211*4b9ad928SBarry Smith EXTERN_C_BEGIN 212*4b9ad928SBarry Smith #undef __FUNCT__ 213*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha_Composite" 214*4b9ad928SBarry Smith int PCCompositeSpecialSetAlpha_Composite(PC pc,PetscScalar alpha) 215*4b9ad928SBarry Smith { 216*4b9ad928SBarry Smith PC_Composite *jac = (PC_Composite*)pc->data; 217*4b9ad928SBarry Smith PetscFunctionBegin; 218*4b9ad928SBarry Smith jac->alpha = alpha; 219*4b9ad928SBarry Smith PetscFunctionReturn(0); 220*4b9ad928SBarry Smith } 221*4b9ad928SBarry Smith EXTERN_C_END 222*4b9ad928SBarry Smith 223*4b9ad928SBarry Smith EXTERN_C_BEGIN 224*4b9ad928SBarry Smith #undef __FUNCT__ 225*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType_Composite" 226*4b9ad928SBarry Smith int PCCompositeSetType_Composite(PC pc,PCCompositeType type) 227*4b9ad928SBarry Smith { 228*4b9ad928SBarry Smith PetscFunctionBegin; 229*4b9ad928SBarry Smith if (type == PC_COMPOSITE_ADDITIVE) { 230*4b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 231*4b9ad928SBarry Smith } else if (type == PC_COMPOSITE_MULTIPLICATIVE) { 232*4b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Multiplicative; 233*4b9ad928SBarry Smith } else if (type == PC_COMPOSITE_SPECIAL) { 234*4b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Special; 235*4b9ad928SBarry Smith } else { 236*4b9ad928SBarry Smith SETERRQ(1,"Unkown composite preconditioner type"); 237*4b9ad928SBarry Smith } 238*4b9ad928SBarry Smith PetscFunctionReturn(0); 239*4b9ad928SBarry Smith } 240*4b9ad928SBarry Smith EXTERN_C_END 241*4b9ad928SBarry Smith 242*4b9ad928SBarry Smith EXTERN_C_BEGIN 243*4b9ad928SBarry Smith #undef __FUNCT__ 244*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC_Composite" 245*4b9ad928SBarry Smith int PCCompositeAddPC_Composite(PC pc,PCType type) 246*4b9ad928SBarry Smith { 247*4b9ad928SBarry Smith PC_Composite *jac; 248*4b9ad928SBarry Smith PC_CompositeLink next,link; 249*4b9ad928SBarry Smith int ierr,cnt = 0; 250*4b9ad928SBarry Smith char *prefix,newprefix[8]; 251*4b9ad928SBarry Smith 252*4b9ad928SBarry Smith PetscFunctionBegin; 253*4b9ad928SBarry Smith ierr = PetscNew(struct _PC_CompositeLink,&link);CHKERRQ(ierr); 254*4b9ad928SBarry Smith link->next = 0; 255*4b9ad928SBarry Smith ierr = PCCreate(pc->comm,&link->pc);CHKERRQ(ierr); 256*4b9ad928SBarry Smith 257*4b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 258*4b9ad928SBarry Smith next = jac->head; 259*4b9ad928SBarry Smith if (!next) { 260*4b9ad928SBarry Smith jac->head = link; 261*4b9ad928SBarry Smith } else { 262*4b9ad928SBarry Smith cnt++; 263*4b9ad928SBarry Smith while (next->next) { 264*4b9ad928SBarry Smith next = next->next; 265*4b9ad928SBarry Smith cnt++; 266*4b9ad928SBarry Smith } 267*4b9ad928SBarry Smith next->next = link; 268*4b9ad928SBarry Smith } 269*4b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 270*4b9ad928SBarry Smith ierr = PCSetOptionsPrefix(link->pc,prefix);CHKERRQ(ierr); 271*4b9ad928SBarry Smith sprintf(newprefix,"sub_%d_",cnt); 272*4b9ad928SBarry Smith ierr = PCAppendOptionsPrefix(link->pc,newprefix);CHKERRQ(ierr); 273*4b9ad928SBarry Smith /* type is set after prefix, because some methods may modify prefix, e.g. pcksp */ 274*4b9ad928SBarry Smith ierr = PCSetType(link->pc,type);CHKERRQ(ierr); 275*4b9ad928SBarry Smith 276*4b9ad928SBarry Smith PetscFunctionReturn(0); 277*4b9ad928SBarry Smith } 278*4b9ad928SBarry Smith EXTERN_C_END 279*4b9ad928SBarry Smith 280*4b9ad928SBarry Smith EXTERN_C_BEGIN 281*4b9ad928SBarry Smith #undef __FUNCT__ 282*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC_Composite" 283*4b9ad928SBarry Smith int PCCompositeGetPC_Composite(PC pc,int n,PC *subpc) 284*4b9ad928SBarry Smith { 285*4b9ad928SBarry Smith PC_Composite *jac; 286*4b9ad928SBarry Smith PC_CompositeLink next; 287*4b9ad928SBarry Smith int i; 288*4b9ad928SBarry Smith 289*4b9ad928SBarry Smith PetscFunctionBegin; 290*4b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 291*4b9ad928SBarry Smith next = jac->head; 292*4b9ad928SBarry Smith for (i=0; i<n; i++) { 293*4b9ad928SBarry Smith if (!next->next) { 294*4b9ad928SBarry Smith SETERRQ(1,"Not enough PCs in composite preconditioner"); 295*4b9ad928SBarry Smith } 296*4b9ad928SBarry Smith next = next->next; 297*4b9ad928SBarry Smith } 298*4b9ad928SBarry Smith *subpc = next->pc; 299*4b9ad928SBarry Smith PetscFunctionReturn(0); 300*4b9ad928SBarry Smith } 301*4b9ad928SBarry Smith EXTERN_C_END 302*4b9ad928SBarry Smith 303*4b9ad928SBarry Smith EXTERN_C_BEGIN 304*4b9ad928SBarry Smith #undef __FUNCT__ 305*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue_Composite" 306*4b9ad928SBarry Smith int PCCompositeSetUseTrue_Composite(PC pc) 307*4b9ad928SBarry Smith { 308*4b9ad928SBarry Smith PC_Composite *jac; 309*4b9ad928SBarry Smith 310*4b9ad928SBarry Smith PetscFunctionBegin; 311*4b9ad928SBarry Smith jac = (PC_Composite*)pc->data; 312*4b9ad928SBarry Smith jac->use_true_matrix = PETSC_TRUE; 313*4b9ad928SBarry Smith PetscFunctionReturn(0); 314*4b9ad928SBarry Smith } 315*4b9ad928SBarry Smith EXTERN_C_END 316*4b9ad928SBarry Smith 317*4b9ad928SBarry Smith /* -------------------------------------------------------------------------------- */ 318*4b9ad928SBarry Smith #undef __FUNCT__ 319*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetType" 320*4b9ad928SBarry Smith /*@C 321*4b9ad928SBarry Smith PCCompositeSetType - Sets the type of composite preconditioner. 322*4b9ad928SBarry Smith 323*4b9ad928SBarry Smith Collective on PC 324*4b9ad928SBarry Smith 325*4b9ad928SBarry Smith Input Parameter: 326*4b9ad928SBarry Smith . pc - the preconditioner context 327*4b9ad928SBarry Smith . type - PC_COMPOSITE_ADDITIVE (default), PC_COMPOSITE_MULTIPLICATIVE, PC_COMPOSITE_SPECIAL 328*4b9ad928SBarry Smith 329*4b9ad928SBarry Smith Options Database Key: 330*4b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 331*4b9ad928SBarry Smith 332*4b9ad928SBarry Smith Level: Developer 333*4b9ad928SBarry Smith 334*4b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 335*4b9ad928SBarry Smith @*/ 336*4b9ad928SBarry Smith int PCCompositeSetType(PC pc,PCCompositeType type) 337*4b9ad928SBarry Smith { 338*4b9ad928SBarry Smith int ierr,(*f)(PC,PCCompositeType); 339*4b9ad928SBarry Smith 340*4b9ad928SBarry Smith PetscFunctionBegin; 341*4b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 342*4b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetType_C",(void (**)(void))&f);CHKERRQ(ierr); 343*4b9ad928SBarry Smith if (f) { 344*4b9ad928SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 345*4b9ad928SBarry Smith } 346*4b9ad928SBarry Smith PetscFunctionReturn(0); 347*4b9ad928SBarry Smith } 348*4b9ad928SBarry Smith 349*4b9ad928SBarry Smith #undef __FUNCT__ 350*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSpecialSetAlpha" 351*4b9ad928SBarry Smith /*@C 352*4b9ad928SBarry Smith PCCompositeSpecialSetAlpha - Sets alpha for the special composite preconditioner 353*4b9ad928SBarry Smith for alphaI + R + S 354*4b9ad928SBarry Smith 355*4b9ad928SBarry Smith Collective on PC 356*4b9ad928SBarry Smith 357*4b9ad928SBarry Smith Input Parameter: 358*4b9ad928SBarry Smith + pc - the preconditioner context 359*4b9ad928SBarry Smith - alpha - scale on identity 360*4b9ad928SBarry Smith 361*4b9ad928SBarry Smith Level: Developer 362*4b9ad928SBarry Smith 363*4b9ad928SBarry Smith .keywords: PC, set, type, composite preconditioner, additive, multiplicative 364*4b9ad928SBarry Smith @*/ 365*4b9ad928SBarry Smith int PCCompositeSpecialSetAlpha(PC pc,PetscScalar alpha) 366*4b9ad928SBarry Smith { 367*4b9ad928SBarry Smith int ierr,(*f)(PC,PetscScalar); 368*4b9ad928SBarry Smith 369*4b9ad928SBarry Smith PetscFunctionBegin; 370*4b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 371*4b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSpecialSetAlpha_C",(void (**)(void))&f);CHKERRQ(ierr); 372*4b9ad928SBarry Smith if (f) { 373*4b9ad928SBarry Smith ierr = (*f)(pc,alpha);CHKERRQ(ierr); 374*4b9ad928SBarry Smith } 375*4b9ad928SBarry Smith PetscFunctionReturn(0); 376*4b9ad928SBarry Smith } 377*4b9ad928SBarry Smith 378*4b9ad928SBarry Smith #undef __FUNCT__ 379*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeAddPC" 380*4b9ad928SBarry Smith /*@C 381*4b9ad928SBarry Smith PCCompositeAddPC - Adds another PC to the composite PC. 382*4b9ad928SBarry Smith 383*4b9ad928SBarry Smith Collective on PC 384*4b9ad928SBarry Smith 385*4b9ad928SBarry Smith Input Parameters: 386*4b9ad928SBarry Smith . pc - the preconditioner context 387*4b9ad928SBarry Smith . type - the type of the new preconditioner 388*4b9ad928SBarry Smith 389*4b9ad928SBarry Smith Level: Developer 390*4b9ad928SBarry Smith 391*4b9ad928SBarry Smith .keywords: PC, composite preconditioner, add 392*4b9ad928SBarry Smith @*/ 393*4b9ad928SBarry Smith int PCCompositeAddPC(PC pc,PCType type) 394*4b9ad928SBarry Smith { 395*4b9ad928SBarry Smith int ierr,(*f)(PC,PCType); 396*4b9ad928SBarry Smith 397*4b9ad928SBarry Smith PetscFunctionBegin; 398*4b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 399*4b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeAddPC_C",(void (**)(void))&f);CHKERRQ(ierr); 400*4b9ad928SBarry Smith if (f) { 401*4b9ad928SBarry Smith ierr = (*f)(pc,type);CHKERRQ(ierr); 402*4b9ad928SBarry Smith } 403*4b9ad928SBarry Smith PetscFunctionReturn(0); 404*4b9ad928SBarry Smith } 405*4b9ad928SBarry Smith 406*4b9ad928SBarry Smith #undef __FUNCT__ 407*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeGetPC" 408*4b9ad928SBarry Smith /*@C 409*4b9ad928SBarry Smith PCCompositeGetPC - Gets one of the PC objects in the composite PC. 410*4b9ad928SBarry Smith 411*4b9ad928SBarry Smith Not Collective 412*4b9ad928SBarry Smith 413*4b9ad928SBarry Smith Input Parameter: 414*4b9ad928SBarry Smith . pc - the preconditioner context 415*4b9ad928SBarry Smith . n - the number of the pc requested 416*4b9ad928SBarry Smith 417*4b9ad928SBarry Smith Output Parameters: 418*4b9ad928SBarry Smith . subpc - the PC requested 419*4b9ad928SBarry Smith 420*4b9ad928SBarry Smith Level: Developer 421*4b9ad928SBarry Smith 422*4b9ad928SBarry Smith .keywords: PC, get, composite preconditioner, sub preconditioner 423*4b9ad928SBarry Smith 424*4b9ad928SBarry Smith .seealso: PCCompositeAddPC() 425*4b9ad928SBarry Smith @*/ 426*4b9ad928SBarry Smith int PCCompositeGetPC(PC pc,int n,PC *subpc) 427*4b9ad928SBarry Smith { 428*4b9ad928SBarry Smith int ierr,(*f)(PC,int,PC *); 429*4b9ad928SBarry Smith 430*4b9ad928SBarry Smith PetscFunctionBegin; 431*4b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 432*4b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 433*4b9ad928SBarry Smith if (f) { 434*4b9ad928SBarry Smith ierr = (*f)(pc,n,subpc);CHKERRQ(ierr); 435*4b9ad928SBarry Smith } else { 436*4b9ad928SBarry Smith SETERRQ(1,"Cannot get pc, not composite type"); 437*4b9ad928SBarry Smith } 438*4b9ad928SBarry Smith PetscFunctionReturn(0); 439*4b9ad928SBarry Smith } 440*4b9ad928SBarry Smith 441*4b9ad928SBarry Smith #undef __FUNCT__ 442*4b9ad928SBarry Smith #define __FUNCT__ "PCCompositeSetUseTrue" 443*4b9ad928SBarry Smith /*@ 444*4b9ad928SBarry Smith PCCompositeSetUseTrue - Sets a flag to indicate that the true matrix (rather than 445*4b9ad928SBarry Smith the matrix used to define the preconditioner) is used to compute 446*4b9ad928SBarry Smith the residual when the multiplicative scheme is used. 447*4b9ad928SBarry Smith 448*4b9ad928SBarry Smith Collective on PC 449*4b9ad928SBarry Smith 450*4b9ad928SBarry Smith Input Parameters: 451*4b9ad928SBarry Smith . pc - the preconditioner context 452*4b9ad928SBarry Smith 453*4b9ad928SBarry Smith Options Database Key: 454*4b9ad928SBarry Smith . -pc_composite_true - Activates PCCompositeSetUseTrue() 455*4b9ad928SBarry Smith 456*4b9ad928SBarry Smith Note: 457*4b9ad928SBarry Smith For the common case in which the preconditioning and linear 458*4b9ad928SBarry Smith system matrices are identical, this routine is unnecessary. 459*4b9ad928SBarry Smith 460*4b9ad928SBarry Smith Level: Developer 461*4b9ad928SBarry Smith 462*4b9ad928SBarry Smith .keywords: PC, composite preconditioner, set, true, flag 463*4b9ad928SBarry Smith 464*4b9ad928SBarry Smith .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal(), PCKSPSetUseTrue() 465*4b9ad928SBarry Smith @*/ 466*4b9ad928SBarry Smith int PCCompositeSetUseTrue(PC pc) 467*4b9ad928SBarry Smith { 468*4b9ad928SBarry Smith int ierr,(*f)(PC); 469*4b9ad928SBarry Smith 470*4b9ad928SBarry Smith PetscFunctionBegin; 471*4b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 472*4b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCCompositeSetUseTrue_C",(void (**)(void))&f);CHKERRQ(ierr); 473*4b9ad928SBarry Smith if (f) { 474*4b9ad928SBarry Smith ierr = (*f)(pc);CHKERRQ(ierr); 475*4b9ad928SBarry Smith } 476*4b9ad928SBarry Smith PetscFunctionReturn(0); 477*4b9ad928SBarry Smith } 478*4b9ad928SBarry Smith 479*4b9ad928SBarry Smith /* -------------------------------------------------------------------------------------------*/ 480*4b9ad928SBarry Smith 481*4b9ad928SBarry Smith /*MC 482*4b9ad928SBarry Smith PCCOMPOSITE - Build a preconditioner by composing together several preconditioners 483*4b9ad928SBarry Smith 484*4b9ad928SBarry Smith Options Database Keys: 485*4b9ad928SBarry Smith . -pc_composite_type <type: one of multiplicative, additive, special> - Sets composite preconditioner type 486*4b9ad928SBarry Smith . -pc_composite_true - Activates PCCompositeSetUseTrue() 487*4b9ad928SBarry Smith 488*4b9ad928SBarry Smith Level: intermediate 489*4b9ad928SBarry Smith 490*4b9ad928SBarry Smith Concepts: composing solvers 491*4b9ad928SBarry Smith 492*4b9ad928SBarry Smith Notes: To use a Krylov method inside the composite preconditioner, set the PCType of one or more 493*4b9ad928SBarry Smith inner PCs to be PCKSP. 494*4b9ad928SBarry Smith Using a Krylov method inside another Krylov method can be dangerous (you get divergence or 495*4b9ad928SBarry Smith the incorrect answer) unless you use KSPFGMRES as the other Krylov method 496*4b9ad928SBarry Smith 497*4b9ad928SBarry Smith 498*4b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 499*4b9ad928SBarry Smith PCSHELL, PCKSP, PCCompositeSetType(), PCCompositeSpecialSetAlpha(), PCCompositeAddPC(), 500*4b9ad928SBarry Smith PCCompositeGetPC(), PCCompositeSetUseTrue() 501*4b9ad928SBarry Smith 502*4b9ad928SBarry Smith M*/ 503*4b9ad928SBarry Smith 504*4b9ad928SBarry Smith EXTERN_C_BEGIN 505*4b9ad928SBarry Smith #undef __FUNCT__ 506*4b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Composite" 507*4b9ad928SBarry Smith int PCCreate_Composite(PC pc) 508*4b9ad928SBarry Smith { 509*4b9ad928SBarry Smith int ierr; 510*4b9ad928SBarry Smith PC_Composite *jac; 511*4b9ad928SBarry Smith 512*4b9ad928SBarry Smith PetscFunctionBegin; 513*4b9ad928SBarry Smith ierr = PetscNew(PC_Composite,&jac);CHKERRQ(ierr); 514*4b9ad928SBarry Smith PetscLogObjectMemory(pc,sizeof(PC_Composite)); 515*4b9ad928SBarry Smith pc->ops->apply = PCApply_Composite_Additive; 516*4b9ad928SBarry Smith pc->ops->setup = PCSetUp_Composite; 517*4b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Composite; 518*4b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Composite; 519*4b9ad928SBarry Smith pc->ops->view = PCView_Composite; 520*4b9ad928SBarry Smith pc->ops->applyrichardson = 0; 521*4b9ad928SBarry Smith 522*4b9ad928SBarry Smith pc->data = (void*)jac; 523*4b9ad928SBarry Smith jac->type = PC_COMPOSITE_ADDITIVE; 524*4b9ad928SBarry Smith jac->work1 = 0; 525*4b9ad928SBarry Smith jac->work2 = 0; 526*4b9ad928SBarry Smith jac->head = 0; 527*4b9ad928SBarry Smith 528*4b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetType_C","PCCompositeSetType_Composite", 529*4b9ad928SBarry Smith PCCompositeSetType_Composite);CHKERRQ(ierr); 530*4b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeAddPC_C","PCCompositeAddPC_Composite", 531*4b9ad928SBarry Smith PCCompositeAddPC_Composite);CHKERRQ(ierr); 532*4b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeGetPC_C","PCCompositeGetPC_Composite", 533*4b9ad928SBarry Smith PCCompositeGetPC_Composite);CHKERRQ(ierr); 534*4b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSetUseTrue_C","PCCompositeSetUseTrue_Composite", 535*4b9ad928SBarry Smith PCCompositeSetUseTrue_Composite);CHKERRQ(ierr); 536*4b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCCompositeSpecialSetAlpha_C","PCCompositeSpecialSetAlpha_Composite", 537*4b9ad928SBarry Smith PCCompositeSpecialSetAlpha_Composite);CHKERRQ(ierr); 538*4b9ad928SBarry Smith 539*4b9ad928SBarry Smith PetscFunctionReturn(0); 540*4b9ad928SBarry Smith } 541*4b9ad928SBarry Smith EXTERN_C_END 542*4b9ad928SBarry Smith 543