1*2a6744ebSBarry Smith #define PETSCKSP_DLL 2*2a6744ebSBarry Smith 3*2a6744ebSBarry Smith /* 4*2a6744ebSBarry Smith Defines a preconditioner defined by R^T S R 5*2a6744ebSBarry Smith */ 6*2a6744ebSBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 7*2a6744ebSBarry Smith #include "petscksp.h" /*I "petscksp.h" I*/ 8*2a6744ebSBarry Smith 9*2a6744ebSBarry Smith typedef struct { 10*2a6744ebSBarry Smith KSP ksp; 11*2a6744ebSBarry Smith Mat R,P; 12*2a6744ebSBarry Smith Vec b,x; 13*2a6744ebSBarry Smith } PC_Galerkin; 14*2a6744ebSBarry Smith 15*2a6744ebSBarry Smith #undef __FUNCT__ 16*2a6744ebSBarry Smith #define __FUNCT__ "PCApply_Galerkin" 17*2a6744ebSBarry Smith static PetscErrorCode PCApply_Galerkin(PC pc,Vec x,Vec y) 18*2a6744ebSBarry Smith { 19*2a6744ebSBarry Smith PetscErrorCode ierr; 20*2a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 21*2a6744ebSBarry Smith 22*2a6744ebSBarry Smith PetscFunctionBegin; 23*2a6744ebSBarry Smith if (jac->R) {ierr = MatRestrict(jac->R,x,jac->b);CHKERRQ(ierr);} 24*2a6744ebSBarry Smith else {ierr = MatRestrict(jac->P,x,jac->b);CHKERRQ(ierr);} 25*2a6744ebSBarry Smith ierr = KSPSolve(jac->ksp,jac->b,jac->x);CHKERRQ(ierr); 26*2a6744ebSBarry Smith if (jac->P) {ierr = MatInterpolate(jac->P,jac->x,y);CHKERRQ(ierr);} 27*2a6744ebSBarry Smith else {ierr = MatInterpolate(jac->R,jac->x,y);CHKERRQ(ierr);} 28*2a6744ebSBarry Smith PetscFunctionReturn(0); 29*2a6744ebSBarry Smith } 30*2a6744ebSBarry Smith 31*2a6744ebSBarry Smith #undef __FUNCT__ 32*2a6744ebSBarry Smith #define __FUNCT__ "PCSetUp_Galerkin" 33*2a6744ebSBarry Smith static PetscErrorCode PCSetUp_Galerkin(PC pc) 34*2a6744ebSBarry Smith { 35*2a6744ebSBarry Smith PetscErrorCode ierr; 36*2a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 37*2a6744ebSBarry Smith Mat a; 38*2a6744ebSBarry Smith 39*2a6744ebSBarry Smith PetscFunctionBegin; 40*2a6744ebSBarry Smith if (!jac->x) { 41*2a6744ebSBarry Smith ierr = KSPGetOperators(jac->ksp,&a,0,0);CHKERRQ(ierr); 42*2a6744ebSBarry Smith if (!a) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set operator of PCGALERKIN KSP with PCGalerkinGetKSP()/KSPSetOperators()"); 43*2a6744ebSBarry Smith ierr = MatGetVecs(a,&jac->x,&jac->b);CHKERRQ(ierr); 44*2a6744ebSBarry Smith } 45*2a6744ebSBarry Smith if (!jac->R && !jac->P) { 46*2a6744ebSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or interpolation of PCGALERKIN with PCGalerkinSetRestriction/Interpolation()"); 47*2a6744ebSBarry Smith } 48*2a6744ebSBarry Smith /* should check here that sizes of R/P match size of a */ 49*2a6744ebSBarry Smith PetscFunctionReturn(0); 50*2a6744ebSBarry Smith } 51*2a6744ebSBarry Smith 52*2a6744ebSBarry Smith #undef __FUNCT__ 53*2a6744ebSBarry Smith #define __FUNCT__ "PCDestroy_Galerkin" 54*2a6744ebSBarry Smith static PetscErrorCode PCDestroy_Galerkin(PC pc) 55*2a6744ebSBarry Smith { 56*2a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 57*2a6744ebSBarry Smith PetscErrorCode ierr; 58*2a6744ebSBarry Smith 59*2a6744ebSBarry Smith PetscFunctionBegin; 60*2a6744ebSBarry Smith if (jac->R) {ierr = MatDestroy(jac->R);CHKERRQ(ierr);} 61*2a6744ebSBarry Smith if (jac->P) {ierr = MatDestroy(jac->P);CHKERRQ(ierr);} 62*2a6744ebSBarry Smith if (jac->x) {ierr = VecDestroy(jac->x);CHKERRQ(ierr);} 63*2a6744ebSBarry Smith if (jac->b) {ierr = VecDestroy(jac->b);CHKERRQ(ierr);} 64*2a6744ebSBarry Smith ierr = KSPDestroy(jac->ksp);CHKERRQ(ierr); 65*2a6744ebSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 66*2a6744ebSBarry Smith PetscFunctionReturn(0); 67*2a6744ebSBarry Smith } 68*2a6744ebSBarry Smith 69*2a6744ebSBarry Smith #undef __FUNCT__ 70*2a6744ebSBarry Smith #define __FUNCT__ "PCView_Galerkin" 71*2a6744ebSBarry Smith static PetscErrorCode PCView_Galerkin(PC pc,PetscViewer viewer) 72*2a6744ebSBarry Smith { 73*2a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 74*2a6744ebSBarry Smith PetscErrorCode ierr; 75*2a6744ebSBarry Smith PetscTruth iascii; 76*2a6744ebSBarry Smith 77*2a6744ebSBarry Smith PetscFunctionBegin; 78*2a6744ebSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 79*2a6744ebSBarry Smith if (iascii) { 80*2a6744ebSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Galerkin PC\n");CHKERRQ(ierr); 81*2a6744ebSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP on Galerkin follow\n");CHKERRQ(ierr); 82*2a6744ebSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 83*2a6744ebSBarry Smith } else { 84*2a6744ebSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCGalerkin",((PetscObject)viewer)->type_name); 85*2a6744ebSBarry Smith } 86*2a6744ebSBarry Smith ierr = KSPView(jac->ksp,viewer);CHKERRQ(ierr); 87*2a6744ebSBarry Smith PetscFunctionReturn(0); 88*2a6744ebSBarry Smith } 89*2a6744ebSBarry Smith 90*2a6744ebSBarry Smith EXTERN_C_BEGIN 91*2a6744ebSBarry Smith #undef __FUNCT__ 92*2a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinGetKSP_Galerkin" 93*2a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP_Galerkin(PC pc,KSP *ksp) 94*2a6744ebSBarry Smith { 95*2a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 96*2a6744ebSBarry Smith 97*2a6744ebSBarry Smith PetscFunctionBegin; 98*2a6744ebSBarry Smith *ksp = jac->ksp; 99*2a6744ebSBarry Smith PetscFunctionReturn(0); 100*2a6744ebSBarry Smith } 101*2a6744ebSBarry Smith EXTERN_C_END 102*2a6744ebSBarry Smith 103*2a6744ebSBarry Smith EXTERN_C_BEGIN 104*2a6744ebSBarry Smith #undef __FUNCT__ 105*2a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinSetRestriction_Galerkin" 106*2a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetRestriction_Galerkin(PC pc,Mat R) 107*2a6744ebSBarry Smith { 108*2a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 109*2a6744ebSBarry Smith PetscErrorCode ierr; 110*2a6744ebSBarry Smith 111*2a6744ebSBarry Smith PetscFunctionBegin; 112*2a6744ebSBarry Smith if (jac->R) {ierr = MatDestroy(jac->R);CHKERRQ(ierr);} 113*2a6744ebSBarry Smith jac->R = R; 114*2a6744ebSBarry Smith ierr = PetscObjectReference((PetscObject)R);CHKERRQ(ierr); 115*2a6744ebSBarry Smith PetscFunctionReturn(0); 116*2a6744ebSBarry Smith } 117*2a6744ebSBarry Smith EXTERN_C_END 118*2a6744ebSBarry Smith 119*2a6744ebSBarry Smith EXTERN_C_BEGIN 120*2a6744ebSBarry Smith #undef __FUNCT__ 121*2a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinSetInterpolation_Galerkin" 122*2a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetInterpolation_Galerkin(PC pc,Mat P) 123*2a6744ebSBarry Smith { 124*2a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 125*2a6744ebSBarry Smith PetscErrorCode ierr; 126*2a6744ebSBarry Smith 127*2a6744ebSBarry Smith PetscFunctionBegin; 128*2a6744ebSBarry Smith if (jac->P) {ierr = MatDestroy(jac->P);CHKERRQ(ierr);} 129*2a6744ebSBarry Smith jac->P = P; 130*2a6744ebSBarry Smith ierr = PetscObjectReference((PetscObject)P);CHKERRQ(ierr); 131*2a6744ebSBarry Smith PetscFunctionReturn(0); 132*2a6744ebSBarry Smith } 133*2a6744ebSBarry Smith EXTERN_C_END 134*2a6744ebSBarry Smith 135*2a6744ebSBarry Smith /* -------------------------------------------------------------------------------- */ 136*2a6744ebSBarry Smith #undef __FUNCT__ 137*2a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinSetRestriction" 138*2a6744ebSBarry Smith /*@ 139*2a6744ebSBarry Smith PCGalerkinSetRestriction - Sets the restriction operator for the "Galerkin-type" preconditioner 140*2a6744ebSBarry Smith 141*2a6744ebSBarry Smith Collective on PC 142*2a6744ebSBarry Smith 143*2a6744ebSBarry Smith Input Parameter: 144*2a6744ebSBarry Smith + pc - the preconditioner context 145*2a6744ebSBarry Smith - R - the restriction operator 146*2a6744ebSBarry Smith 147*2a6744ebSBarry Smith Notes: Either this or PCGalerkinSetInterpolation() or both must be called 148*2a6744ebSBarry Smith 149*2a6744ebSBarry Smith Level: Intermediate 150*2a6744ebSBarry Smith 151*2a6744ebSBarry Smith .keywords: PC, set, Galerkin preconditioner 152*2a6744ebSBarry Smith 153*2a6744ebSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 154*2a6744ebSBarry Smith PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 155*2a6744ebSBarry Smith 156*2a6744ebSBarry Smith @*/ 157*2a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetRestriction(PC pc,Mat R) 158*2a6744ebSBarry Smith { 159*2a6744ebSBarry Smith PetscErrorCode ierr,(*f)(PC,Mat); 160*2a6744ebSBarry Smith 161*2a6744ebSBarry Smith PetscFunctionBegin; 162*2a6744ebSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 163*2a6744ebSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCGalerkinSetRestriction_C",(void (**)(void))&f);CHKERRQ(ierr); 164*2a6744ebSBarry Smith if (f) { 165*2a6744ebSBarry Smith ierr = (*f)(pc,R);CHKERRQ(ierr); 166*2a6744ebSBarry Smith } 167*2a6744ebSBarry Smith PetscFunctionReturn(0); 168*2a6744ebSBarry Smith } 169*2a6744ebSBarry Smith 170*2a6744ebSBarry Smith #undef __FUNCT__ 171*2a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinSetRestriction" 172*2a6744ebSBarry Smith /*@ 173*2a6744ebSBarry Smith PCGalerkinSetInterpolation - Sets the interpolation operator for the "Galerkin-type" preconditioner 174*2a6744ebSBarry Smith 175*2a6744ebSBarry Smith Collective on PC 176*2a6744ebSBarry Smith 177*2a6744ebSBarry Smith Input Parameter: 178*2a6744ebSBarry Smith + pc - the preconditioner context 179*2a6744ebSBarry Smith - R - the interpolation operator 180*2a6744ebSBarry Smith 181*2a6744ebSBarry Smith Notes: Either this or PCGalerkinSetRestriction() or both must be called 182*2a6744ebSBarry Smith 183*2a6744ebSBarry Smith Level: Intermediate 184*2a6744ebSBarry Smith 185*2a6744ebSBarry Smith .keywords: PC, set, Galerkin preconditioner 186*2a6744ebSBarry Smith 187*2a6744ebSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 188*2a6744ebSBarry Smith PCGalerkinSetRestriction(), PCGalerkinGetKSP() 189*2a6744ebSBarry Smith 190*2a6744ebSBarry Smith @*/ 191*2a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetInterpolation(PC pc,Mat P) 192*2a6744ebSBarry Smith { 193*2a6744ebSBarry Smith PetscErrorCode ierr,(*f)(PC,Mat); 194*2a6744ebSBarry Smith 195*2a6744ebSBarry Smith PetscFunctionBegin; 196*2a6744ebSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 197*2a6744ebSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCGalerkinSetInterpolation_C",(void (**)(void))&f);CHKERRQ(ierr); 198*2a6744ebSBarry Smith if (f) { 199*2a6744ebSBarry Smith ierr = (*f)(pc,P);CHKERRQ(ierr); 200*2a6744ebSBarry Smith } 201*2a6744ebSBarry Smith PetscFunctionReturn(0); 202*2a6744ebSBarry Smith } 203*2a6744ebSBarry Smith 204*2a6744ebSBarry Smith #undef __FUNCT__ 205*2a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinGetKSP" 206*2a6744ebSBarry Smith /*@ 207*2a6744ebSBarry Smith PCGalerkinGetKSP - Gets the KSP object in the Galerkin PC. 208*2a6744ebSBarry Smith 209*2a6744ebSBarry Smith Not Collective 210*2a6744ebSBarry Smith 211*2a6744ebSBarry Smith Input Parameter: 212*2a6744ebSBarry Smith . pc - the preconditioner context 213*2a6744ebSBarry Smith 214*2a6744ebSBarry Smith Output Parameters: 215*2a6744ebSBarry Smith . ksp - the KSP object 216*2a6744ebSBarry Smith 217*2a6744ebSBarry Smith Level: Intermediate 218*2a6744ebSBarry Smith 219*2a6744ebSBarry Smith .keywords: PC, get, Galerkin preconditioner, sub preconditioner 220*2a6744ebSBarry Smith 221*2a6744ebSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 222*2a6744ebSBarry Smith PCGalerkinSetRestriction(), PCGalerkinSetInterpolation() 223*2a6744ebSBarry Smith 224*2a6744ebSBarry Smith @*/ 225*2a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC pc,KSP *ksp) 226*2a6744ebSBarry Smith { 227*2a6744ebSBarry Smith PetscErrorCode ierr,(*f)(PC,KSP *); 228*2a6744ebSBarry Smith 229*2a6744ebSBarry Smith PetscFunctionBegin; 230*2a6744ebSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 231*2a6744ebSBarry Smith PetscValidPointer(ksp,2); 232*2a6744ebSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCGalerkinGetKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 233*2a6744ebSBarry Smith if (f) { 234*2a6744ebSBarry Smith ierr = (*f)(pc,ksp);CHKERRQ(ierr); 235*2a6744ebSBarry Smith } else { 236*2a6744ebSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get KSP, not Galerkin type"); 237*2a6744ebSBarry Smith } 238*2a6744ebSBarry Smith PetscFunctionReturn(0); 239*2a6744ebSBarry Smith } 240*2a6744ebSBarry Smith 241*2a6744ebSBarry Smith 242*2a6744ebSBarry Smith /* -------------------------------------------------------------------------------------------*/ 243*2a6744ebSBarry Smith 244*2a6744ebSBarry Smith /*MC 245*2a6744ebSBarry Smith PCGALERKIN - Build (part of) a preconditioner by P S R (where P is often R^T) 246*2a6744ebSBarry Smith 247*2a6744ebSBarry Smith $ Use PCGalerkinSetRestriction(pc,R) and/or PCGalerkinSetInterpolation(pc,P) followed by 248*2a6744ebSBarry Smith $ PCGalerkinGetKSP(pc,&ksp); KSPSetOperations(ksp,A,....) 249*2a6744ebSBarry Smith 250*2a6744ebSBarry Smith Level: intermediate 251*2a6744ebSBarry Smith 252*2a6744ebSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 253*2a6744ebSBarry Smith PCSHELL, PCKSP, PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 254*2a6744ebSBarry Smith 255*2a6744ebSBarry Smith M*/ 256*2a6744ebSBarry Smith 257*2a6744ebSBarry Smith EXTERN_C_BEGIN 258*2a6744ebSBarry Smith #undef __FUNCT__ 259*2a6744ebSBarry Smith #define __FUNCT__ "PCCreate_Galerkin" 260*2a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Galerkin(PC pc) 261*2a6744ebSBarry Smith { 262*2a6744ebSBarry Smith PetscErrorCode ierr; 263*2a6744ebSBarry Smith PC_Galerkin *jac; 264*2a6744ebSBarry Smith 265*2a6744ebSBarry Smith PetscFunctionBegin; 266*2a6744ebSBarry Smith ierr = PetscNew(PC_Galerkin,&jac);CHKERRQ(ierr); 267*2a6744ebSBarry Smith pc->ops->apply = PCApply_Galerkin; 268*2a6744ebSBarry Smith pc->ops->setup = PCSetUp_Galerkin; 269*2a6744ebSBarry Smith pc->ops->destroy = PCDestroy_Galerkin; 270*2a6744ebSBarry Smith pc->ops->view = PCView_Galerkin; 271*2a6744ebSBarry Smith pc->ops->applyrichardson = 0; 272*2a6744ebSBarry Smith 273*2a6744ebSBarry Smith ierr = KSPCreate(pc->comm,&jac->ksp);CHKERRQ(ierr); 274*2a6744ebSBarry Smith 275*2a6744ebSBarry Smith pc->data = (void*)jac; 276*2a6744ebSBarry Smith 277*2a6744ebSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGalerkinSetRestriction_C","PCGalerkinSetRestriction_Galerkin", 278*2a6744ebSBarry Smith PCGalerkinSetRestriction_Galerkin);CHKERRQ(ierr); 279*2a6744ebSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGalerkinSetInterpolation_C","PCGalerkinSetInterpolation_Galerkin", 280*2a6744ebSBarry Smith PCGalerkinSetInterpolation_Galerkin);CHKERRQ(ierr); 281*2a6744ebSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGalerkinGetKSP_C","PCGalerkinGetKSP_Galerkin", 282*2a6744ebSBarry Smith PCGalerkinGetKSP_Galerkin);CHKERRQ(ierr); 283*2a6744ebSBarry Smith PetscFunctionReturn(0); 284*2a6744ebSBarry Smith } 285*2a6744ebSBarry Smith EXTERN_C_END 286*2a6744ebSBarry Smith 287