12a6744ebSBarry Smith #define PETSCKSP_DLL 22a6744ebSBarry Smith 32a6744ebSBarry Smith /* 42a6744ebSBarry Smith Defines a preconditioner defined by R^T S R 52a6744ebSBarry Smith */ 62a6744ebSBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 72a6744ebSBarry Smith #include "petscksp.h" /*I "petscksp.h" I*/ 82a6744ebSBarry Smith 92a6744ebSBarry Smith typedef struct { 102a6744ebSBarry Smith KSP ksp; 112a6744ebSBarry Smith Mat R,P; 122a6744ebSBarry Smith Vec b,x; 132a6744ebSBarry Smith } PC_Galerkin; 142a6744ebSBarry Smith 152a6744ebSBarry Smith #undef __FUNCT__ 162a6744ebSBarry Smith #define __FUNCT__ "PCApply_Galerkin" 172a6744ebSBarry Smith static PetscErrorCode PCApply_Galerkin(PC pc,Vec x,Vec y) 182a6744ebSBarry Smith { 192a6744ebSBarry Smith PetscErrorCode ierr; 202a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 212a6744ebSBarry Smith 222a6744ebSBarry Smith PetscFunctionBegin; 232a6744ebSBarry Smith if (jac->R) {ierr = MatRestrict(jac->R,x,jac->b);CHKERRQ(ierr);} 242a6744ebSBarry Smith else {ierr = MatRestrict(jac->P,x,jac->b);CHKERRQ(ierr);} 252a6744ebSBarry Smith ierr = KSPSolve(jac->ksp,jac->b,jac->x);CHKERRQ(ierr); 262a6744ebSBarry Smith if (jac->P) {ierr = MatInterpolate(jac->P,jac->x,y);CHKERRQ(ierr);} 272a6744ebSBarry Smith else {ierr = MatInterpolate(jac->R,jac->x,y);CHKERRQ(ierr);} 282a6744ebSBarry Smith PetscFunctionReturn(0); 292a6744ebSBarry Smith } 302a6744ebSBarry Smith 312a6744ebSBarry Smith #undef __FUNCT__ 322a6744ebSBarry Smith #define __FUNCT__ "PCSetUp_Galerkin" 332a6744ebSBarry Smith static PetscErrorCode PCSetUp_Galerkin(PC pc) 342a6744ebSBarry Smith { 352a6744ebSBarry Smith PetscErrorCode ierr; 362a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 37*906ed7ccSBarry Smith PetscTruth a; 38*906ed7ccSBarry Smith Vec *xx,*yy; 392a6744ebSBarry Smith 402a6744ebSBarry Smith PetscFunctionBegin; 412a6744ebSBarry Smith if (!jac->x) { 42*906ed7ccSBarry Smith ierr = KSPGetOperatorsSet(jac->ksp,&a,PETSC_NULL);CHKERRQ(ierr); 432a6744ebSBarry Smith if (!a) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set operator of PCGALERKIN KSP with PCGalerkinGetKSP()/KSPSetOperators()"); 44*906ed7ccSBarry Smith ierr = KSPGetVecs(jac->ksp,1,&xx,1,&yy);CHKERRQ(ierr); 45*906ed7ccSBarry Smith jac->x = *xx; 46*906ed7ccSBarry Smith jac->b = *yy; 47*906ed7ccSBarry Smith ierr = PetscFree(xx);CHKERRQ(ierr); 48*906ed7ccSBarry Smith ierr = PetscFree(yy);CHKERRQ(ierr); 492a6744ebSBarry Smith } 502a6744ebSBarry Smith if (!jac->R && !jac->P) { 512a6744ebSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or interpolation of PCGALERKIN with PCGalerkinSetRestriction/Interpolation()"); 522a6744ebSBarry Smith } 532a6744ebSBarry Smith /* should check here that sizes of R/P match size of a */ 542a6744ebSBarry Smith PetscFunctionReturn(0); 552a6744ebSBarry Smith } 562a6744ebSBarry Smith 572a6744ebSBarry Smith #undef __FUNCT__ 582a6744ebSBarry Smith #define __FUNCT__ "PCDestroy_Galerkin" 592a6744ebSBarry Smith static PetscErrorCode PCDestroy_Galerkin(PC pc) 602a6744ebSBarry Smith { 612a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 622a6744ebSBarry Smith PetscErrorCode ierr; 632a6744ebSBarry Smith 642a6744ebSBarry Smith PetscFunctionBegin; 652a6744ebSBarry Smith if (jac->R) {ierr = MatDestroy(jac->R);CHKERRQ(ierr);} 662a6744ebSBarry Smith if (jac->P) {ierr = MatDestroy(jac->P);CHKERRQ(ierr);} 672a6744ebSBarry Smith if (jac->x) {ierr = VecDestroy(jac->x);CHKERRQ(ierr);} 682a6744ebSBarry Smith if (jac->b) {ierr = VecDestroy(jac->b);CHKERRQ(ierr);} 692a6744ebSBarry Smith ierr = KSPDestroy(jac->ksp);CHKERRQ(ierr); 702a6744ebSBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 712a6744ebSBarry Smith PetscFunctionReturn(0); 722a6744ebSBarry Smith } 732a6744ebSBarry Smith 742a6744ebSBarry Smith #undef __FUNCT__ 752a6744ebSBarry Smith #define __FUNCT__ "PCView_Galerkin" 762a6744ebSBarry Smith static PetscErrorCode PCView_Galerkin(PC pc,PetscViewer viewer) 772a6744ebSBarry Smith { 782a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 792a6744ebSBarry Smith PetscErrorCode ierr; 802a6744ebSBarry Smith PetscTruth iascii; 812a6744ebSBarry Smith 822a6744ebSBarry Smith PetscFunctionBegin; 832a6744ebSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 842a6744ebSBarry Smith if (iascii) { 852a6744ebSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Galerkin PC\n");CHKERRQ(ierr); 862a6744ebSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"KSP on Galerkin follow\n");CHKERRQ(ierr); 872a6744ebSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 882a6744ebSBarry Smith } else { 892a6744ebSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCGalerkin",((PetscObject)viewer)->type_name); 902a6744ebSBarry Smith } 912a6744ebSBarry Smith ierr = KSPView(jac->ksp,viewer);CHKERRQ(ierr); 922a6744ebSBarry Smith PetscFunctionReturn(0); 932a6744ebSBarry Smith } 942a6744ebSBarry Smith 952a6744ebSBarry Smith EXTERN_C_BEGIN 962a6744ebSBarry Smith #undef __FUNCT__ 972a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinGetKSP_Galerkin" 982a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP_Galerkin(PC pc,KSP *ksp) 992a6744ebSBarry Smith { 1002a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 1012a6744ebSBarry Smith 1022a6744ebSBarry Smith PetscFunctionBegin; 1032a6744ebSBarry Smith *ksp = jac->ksp; 1042a6744ebSBarry Smith PetscFunctionReturn(0); 1052a6744ebSBarry Smith } 1062a6744ebSBarry Smith EXTERN_C_END 1072a6744ebSBarry Smith 1082a6744ebSBarry Smith EXTERN_C_BEGIN 1092a6744ebSBarry Smith #undef __FUNCT__ 1102a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinSetRestriction_Galerkin" 1112a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetRestriction_Galerkin(PC pc,Mat R) 1122a6744ebSBarry Smith { 1132a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 1142a6744ebSBarry Smith PetscErrorCode ierr; 1152a6744ebSBarry Smith 1162a6744ebSBarry Smith PetscFunctionBegin; 1172a6744ebSBarry Smith if (jac->R) {ierr = MatDestroy(jac->R);CHKERRQ(ierr);} 1182a6744ebSBarry Smith jac->R = R; 1192a6744ebSBarry Smith ierr = PetscObjectReference((PetscObject)R);CHKERRQ(ierr); 1202a6744ebSBarry Smith PetscFunctionReturn(0); 1212a6744ebSBarry Smith } 1222a6744ebSBarry Smith EXTERN_C_END 1232a6744ebSBarry Smith 1242a6744ebSBarry Smith EXTERN_C_BEGIN 1252a6744ebSBarry Smith #undef __FUNCT__ 1262a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinSetInterpolation_Galerkin" 1272a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetInterpolation_Galerkin(PC pc,Mat P) 1282a6744ebSBarry Smith { 1292a6744ebSBarry Smith PC_Galerkin *jac = (PC_Galerkin*)pc->data; 1302a6744ebSBarry Smith PetscErrorCode ierr; 1312a6744ebSBarry Smith 1322a6744ebSBarry Smith PetscFunctionBegin; 1332a6744ebSBarry Smith if (jac->P) {ierr = MatDestroy(jac->P);CHKERRQ(ierr);} 1342a6744ebSBarry Smith jac->P = P; 1352a6744ebSBarry Smith ierr = PetscObjectReference((PetscObject)P);CHKERRQ(ierr); 1362a6744ebSBarry Smith PetscFunctionReturn(0); 1372a6744ebSBarry Smith } 1382a6744ebSBarry Smith EXTERN_C_END 1392a6744ebSBarry Smith 1402a6744ebSBarry Smith /* -------------------------------------------------------------------------------- */ 1412a6744ebSBarry Smith #undef __FUNCT__ 1422a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinSetRestriction" 1432a6744ebSBarry Smith /*@ 1442a6744ebSBarry Smith PCGalerkinSetRestriction - Sets the restriction operator for the "Galerkin-type" preconditioner 1452a6744ebSBarry Smith 1462a6744ebSBarry Smith Collective on PC 1472a6744ebSBarry Smith 1482a6744ebSBarry Smith Input Parameter: 1492a6744ebSBarry Smith + pc - the preconditioner context 1502a6744ebSBarry Smith - R - the restriction operator 1512a6744ebSBarry Smith 1522a6744ebSBarry Smith Notes: Either this or PCGalerkinSetInterpolation() or both must be called 1532a6744ebSBarry Smith 1542a6744ebSBarry Smith Level: Intermediate 1552a6744ebSBarry Smith 1562a6744ebSBarry Smith .keywords: PC, set, Galerkin preconditioner 1572a6744ebSBarry Smith 1582a6744ebSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 1592a6744ebSBarry Smith PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 1602a6744ebSBarry Smith 1612a6744ebSBarry Smith @*/ 1622a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetRestriction(PC pc,Mat R) 1632a6744ebSBarry Smith { 1642a6744ebSBarry Smith PetscErrorCode ierr,(*f)(PC,Mat); 1652a6744ebSBarry Smith 1662a6744ebSBarry Smith PetscFunctionBegin; 1672a6744ebSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 1682a6744ebSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCGalerkinSetRestriction_C",(void (**)(void))&f);CHKERRQ(ierr); 1692a6744ebSBarry Smith if (f) { 1702a6744ebSBarry Smith ierr = (*f)(pc,R);CHKERRQ(ierr); 1712a6744ebSBarry Smith } 1722a6744ebSBarry Smith PetscFunctionReturn(0); 1732a6744ebSBarry Smith } 1742a6744ebSBarry Smith 1752a6744ebSBarry Smith #undef __FUNCT__ 1762a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinSetRestriction" 1772a6744ebSBarry Smith /*@ 1782a6744ebSBarry Smith PCGalerkinSetInterpolation - Sets the interpolation operator for the "Galerkin-type" preconditioner 1792a6744ebSBarry Smith 1802a6744ebSBarry Smith Collective on PC 1812a6744ebSBarry Smith 1822a6744ebSBarry Smith Input Parameter: 1832a6744ebSBarry Smith + pc - the preconditioner context 1842a6744ebSBarry Smith - R - the interpolation operator 1852a6744ebSBarry Smith 1862a6744ebSBarry Smith Notes: Either this or PCGalerkinSetRestriction() or both must be called 1872a6744ebSBarry Smith 1882a6744ebSBarry Smith Level: Intermediate 1892a6744ebSBarry Smith 1902a6744ebSBarry Smith .keywords: PC, set, Galerkin preconditioner 1912a6744ebSBarry Smith 1922a6744ebSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 1932a6744ebSBarry Smith PCGalerkinSetRestriction(), PCGalerkinGetKSP() 1942a6744ebSBarry Smith 1952a6744ebSBarry Smith @*/ 1962a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinSetInterpolation(PC pc,Mat P) 1972a6744ebSBarry Smith { 1982a6744ebSBarry Smith PetscErrorCode ierr,(*f)(PC,Mat); 1992a6744ebSBarry Smith 2002a6744ebSBarry Smith PetscFunctionBegin; 2012a6744ebSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2022a6744ebSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCGalerkinSetInterpolation_C",(void (**)(void))&f);CHKERRQ(ierr); 2032a6744ebSBarry Smith if (f) { 2042a6744ebSBarry Smith ierr = (*f)(pc,P);CHKERRQ(ierr); 2052a6744ebSBarry Smith } 2062a6744ebSBarry Smith PetscFunctionReturn(0); 2072a6744ebSBarry Smith } 2082a6744ebSBarry Smith 2092a6744ebSBarry Smith #undef __FUNCT__ 2102a6744ebSBarry Smith #define __FUNCT__ "PCGalerkinGetKSP" 2112a6744ebSBarry Smith /*@ 2122a6744ebSBarry Smith PCGalerkinGetKSP - Gets the KSP object in the Galerkin PC. 2132a6744ebSBarry Smith 2142a6744ebSBarry Smith Not Collective 2152a6744ebSBarry Smith 2162a6744ebSBarry Smith Input Parameter: 2172a6744ebSBarry Smith . pc - the preconditioner context 2182a6744ebSBarry Smith 2192a6744ebSBarry Smith Output Parameters: 2202a6744ebSBarry Smith . ksp - the KSP object 2212a6744ebSBarry Smith 2222a6744ebSBarry Smith Level: Intermediate 2232a6744ebSBarry Smith 2242a6744ebSBarry Smith .keywords: PC, get, Galerkin preconditioner, sub preconditioner 2252a6744ebSBarry Smith 2262a6744ebSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 2272a6744ebSBarry Smith PCGalerkinSetRestriction(), PCGalerkinSetInterpolation() 2282a6744ebSBarry Smith 2292a6744ebSBarry Smith @*/ 2302a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCGalerkinGetKSP(PC pc,KSP *ksp) 2312a6744ebSBarry Smith { 2322a6744ebSBarry Smith PetscErrorCode ierr,(*f)(PC,KSP *); 2332a6744ebSBarry Smith 2342a6744ebSBarry Smith PetscFunctionBegin; 2352a6744ebSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2362a6744ebSBarry Smith PetscValidPointer(ksp,2); 2372a6744ebSBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCGalerkinGetKSP_C",(void (**)(void))&f);CHKERRQ(ierr); 2382a6744ebSBarry Smith if (f) { 2392a6744ebSBarry Smith ierr = (*f)(pc,ksp);CHKERRQ(ierr); 2402a6744ebSBarry Smith } else { 2412a6744ebSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,"Cannot get KSP, not Galerkin type"); 2422a6744ebSBarry Smith } 2432a6744ebSBarry Smith PetscFunctionReturn(0); 2442a6744ebSBarry Smith } 2452a6744ebSBarry Smith 2462a6744ebSBarry Smith 2472a6744ebSBarry Smith /* -------------------------------------------------------------------------------------------*/ 2482a6744ebSBarry Smith 2492a6744ebSBarry Smith /*MC 2502a6744ebSBarry Smith PCGALERKIN - Build (part of) a preconditioner by P S R (where P is often R^T) 2512a6744ebSBarry Smith 2522a6744ebSBarry Smith $ Use PCGalerkinSetRestriction(pc,R) and/or PCGalerkinSetInterpolation(pc,P) followed by 2532a6744ebSBarry Smith $ PCGalerkinGetKSP(pc,&ksp); KSPSetOperations(ksp,A,....) 2542a6744ebSBarry Smith 2552a6744ebSBarry Smith Level: intermediate 2562a6744ebSBarry Smith 2572a6744ebSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 2582a6744ebSBarry Smith PCSHELL, PCKSP, PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 2592a6744ebSBarry Smith 2602a6744ebSBarry Smith M*/ 2612a6744ebSBarry Smith 2622a6744ebSBarry Smith EXTERN_C_BEGIN 2632a6744ebSBarry Smith #undef __FUNCT__ 2642a6744ebSBarry Smith #define __FUNCT__ "PCCreate_Galerkin" 2652a6744ebSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Galerkin(PC pc) 2662a6744ebSBarry Smith { 2672a6744ebSBarry Smith PetscErrorCode ierr; 2682a6744ebSBarry Smith PC_Galerkin *jac; 2692a6744ebSBarry Smith 2702a6744ebSBarry Smith PetscFunctionBegin; 2712a6744ebSBarry Smith ierr = PetscNew(PC_Galerkin,&jac);CHKERRQ(ierr); 2722a6744ebSBarry Smith pc->ops->apply = PCApply_Galerkin; 2732a6744ebSBarry Smith pc->ops->setup = PCSetUp_Galerkin; 2742a6744ebSBarry Smith pc->ops->destroy = PCDestroy_Galerkin; 2752a6744ebSBarry Smith pc->ops->view = PCView_Galerkin; 2762a6744ebSBarry Smith pc->ops->applyrichardson = 0; 2772a6744ebSBarry Smith 2782a6744ebSBarry Smith ierr = KSPCreate(pc->comm,&jac->ksp);CHKERRQ(ierr); 2792a6744ebSBarry Smith 2802a6744ebSBarry Smith pc->data = (void*)jac; 2812a6744ebSBarry Smith 2822a6744ebSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGalerkinSetRestriction_C","PCGalerkinSetRestriction_Galerkin", 2832a6744ebSBarry Smith PCGalerkinSetRestriction_Galerkin);CHKERRQ(ierr); 2842a6744ebSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGalerkinSetInterpolation_C","PCGalerkinSetInterpolation_Galerkin", 2852a6744ebSBarry Smith PCGalerkinSetInterpolation_Galerkin);CHKERRQ(ierr); 2862a6744ebSBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGalerkinGetKSP_C","PCGalerkinGetKSP_Galerkin", 2872a6744ebSBarry Smith PCGalerkinGetKSP_Galerkin);CHKERRQ(ierr); 2882a6744ebSBarry Smith PetscFunctionReturn(0); 2892a6744ebSBarry Smith } 2902a6744ebSBarry Smith EXTERN_C_END 2912a6744ebSBarry Smith 292