14b9ad928SBarry Smith /*$Id: sor.c,v 1.103 2001/08/21 21:03:17 bsmith Exp $*/ 24b9ad928SBarry Smith 34b9ad928SBarry Smith /* 44b9ad928SBarry Smith Defines a (S)SOR preconditioner for any Mat implementation 54b9ad928SBarry Smith */ 64b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 74b9ad928SBarry Smith 84b9ad928SBarry Smith typedef struct { 94b9ad928SBarry Smith int its; /* inner iterations, number of sweeps */ 104b9ad928SBarry Smith int lits; /* local inner iterations, number of sweeps applied by the local matrix mat->A */ 114b9ad928SBarry Smith MatSORType sym; /* forward, reverse, symmetric etc. */ 124b9ad928SBarry Smith PetscReal omega; 134b9ad928SBarry Smith } PC_SOR; 144b9ad928SBarry Smith 154b9ad928SBarry Smith #undef __FUNCT__ 164b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_SOR" 174b9ad928SBarry Smith static int PCDestroy_SOR(PC pc) 184b9ad928SBarry Smith { 194b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 204b9ad928SBarry Smith int ierr; 214b9ad928SBarry Smith 224b9ad928SBarry Smith PetscFunctionBegin; 234b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 244b9ad928SBarry Smith PetscFunctionReturn(0); 254b9ad928SBarry Smith } 264b9ad928SBarry Smith 274b9ad928SBarry Smith #undef __FUNCT__ 284b9ad928SBarry Smith #define __FUNCT__ "PCApply_SOR" 294b9ad928SBarry Smith static int PCApply_SOR(PC pc,Vec x,Vec y) 304b9ad928SBarry Smith { 314b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 324b9ad928SBarry Smith int ierr,flag = jac->sym | SOR_ZERO_INITIAL_GUESS; 334b9ad928SBarry Smith 344b9ad928SBarry Smith PetscFunctionBegin; 354b9ad928SBarry Smith ierr = MatRelax(pc->pmat,x,jac->omega,(MatSORType)flag,0.0,jac->its,jac->lits,y);CHKERRQ(ierr); 364b9ad928SBarry Smith PetscFunctionReturn(0); 374b9ad928SBarry Smith } 384b9ad928SBarry Smith 394b9ad928SBarry Smith #undef __FUNCT__ 404b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_SOR" 414b9ad928SBarry Smith static int PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal atol, PetscReal dtol,int its) 424b9ad928SBarry Smith { 434b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 444b9ad928SBarry Smith int ierr; 454b9ad928SBarry Smith 464b9ad928SBarry Smith PetscFunctionBegin; 474b9ad928SBarry Smith PetscLogInfo(pc,"PCApplyRichardson_SOR: Warning, convergence critera ignored, using %d iterations\n",its); 484b9ad928SBarry Smith its = its*jac->its; 494b9ad928SBarry Smith ierr = MatRelax(pc->pmat,b,jac->omega,(MatSORType)jac->sym,0.0,its,jac->lits,y);CHKERRQ(ierr); 504b9ad928SBarry Smith PetscFunctionReturn(0); 514b9ad928SBarry Smith } 524b9ad928SBarry Smith 534b9ad928SBarry Smith #undef __FUNCT__ 544b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR" 554b9ad928SBarry Smith static int PCSetFromOptions_SOR(PC pc) 564b9ad928SBarry Smith { 574b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 584b9ad928SBarry Smith int ierr; 594b9ad928SBarry Smith PetscTruth flg; 604b9ad928SBarry Smith 614b9ad928SBarry Smith PetscFunctionBegin; 624b9ad928SBarry Smith ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr); 634b9ad928SBarry Smith ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);CHKERRQ(ierr); 644b9ad928SBarry Smith ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr); 654b9ad928SBarry Smith ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);CHKERRQ(ierr); 664b9ad928SBarry Smith ierr = PetscOptionsLogicalGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 674b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 684b9ad928SBarry Smith ierr = PetscOptionsLogicalGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 694b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);} 704b9ad928SBarry Smith ierr = PetscOptionsLogicalGroup("-pc_sor_local_symmetric","use SSOR seperately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 714b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 724b9ad928SBarry Smith ierr = PetscOptionsLogicalGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 734b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);} 744b9ad928SBarry Smith ierr = PetscOptionsLogicalGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 754b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);} 764b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 774b9ad928SBarry Smith PetscFunctionReturn(0); 784b9ad928SBarry Smith } 794b9ad928SBarry Smith 804b9ad928SBarry Smith #undef __FUNCT__ 814b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR" 824b9ad928SBarry Smith static int PCView_SOR(PC pc,PetscViewer viewer) 834b9ad928SBarry Smith { 844b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 854b9ad928SBarry Smith MatSORType sym = jac->sym; 86*2fc52814SBarry Smith const char *sortype; 874b9ad928SBarry Smith int ierr; 884b9ad928SBarry Smith PetscTruth isascii; 894b9ad928SBarry Smith 904b9ad928SBarry Smith PetscFunctionBegin; 914b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 924b9ad928SBarry Smith if (isascii) { 934b9ad928SBarry Smith if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer," SOR: zero initial guess\n");CHKERRQ(ierr);} 944b9ad928SBarry Smith if (sym == SOR_APPLY_UPPER) sortype = "apply_upper"; 954b9ad928SBarry Smith else if (sym == SOR_APPLY_LOWER) sortype = "apply_lower"; 964b9ad928SBarry Smith else if (sym & SOR_EISENSTAT) sortype = "Eisenstat"; 974b9ad928SBarry Smith else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP) 984b9ad928SBarry Smith sortype = "symmetric"; 994b9ad928SBarry Smith else if (sym & SOR_BACKWARD_SWEEP) sortype = "backward"; 1004b9ad928SBarry Smith else if (sym & SOR_FORWARD_SWEEP) sortype = "forward"; 1014b9ad928SBarry Smith else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) 1024b9ad928SBarry Smith sortype = "local_symmetric"; 1034b9ad928SBarry Smith else if (sym & SOR_LOCAL_FORWARD_SWEEP) sortype = "local_forward"; 1044b9ad928SBarry Smith else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward"; 1054b9ad928SBarry Smith else sortype = "unknown"; 1064b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SOR: type = %s, iterations = %d, omega = %g\n",sortype,jac->its,jac->omega);CHKERRQ(ierr); 1074b9ad928SBarry Smith } else { 1084b9ad928SBarry Smith SETERRQ1(1,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name); 1094b9ad928SBarry Smith } 1104b9ad928SBarry Smith PetscFunctionReturn(0); 1114b9ad928SBarry Smith } 1124b9ad928SBarry Smith 1134b9ad928SBarry Smith 1144b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1154b9ad928SBarry Smith EXTERN_C_BEGIN 1164b9ad928SBarry Smith #undef __FUNCT__ 1174b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR" 1184b9ad928SBarry Smith int PCSORSetSymmetric_SOR(PC pc,MatSORType flag) 1194b9ad928SBarry Smith { 1204b9ad928SBarry Smith PC_SOR *jac; 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith PetscFunctionBegin; 1234b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 1244b9ad928SBarry Smith jac->sym = flag; 1254b9ad928SBarry Smith PetscFunctionReturn(0); 1264b9ad928SBarry Smith } 1274b9ad928SBarry Smith EXTERN_C_END 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith EXTERN_C_BEGIN 1304b9ad928SBarry Smith #undef __FUNCT__ 1314b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR" 1324b9ad928SBarry Smith int PCSORSetOmega_SOR(PC pc,PetscReal omega) 1334b9ad928SBarry Smith { 1344b9ad928SBarry Smith PC_SOR *jac; 1354b9ad928SBarry Smith 1364b9ad928SBarry Smith PetscFunctionBegin; 1374b9ad928SBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 1384b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 1394b9ad928SBarry Smith jac->omega = omega; 1404b9ad928SBarry Smith PetscFunctionReturn(0); 1414b9ad928SBarry Smith } 1424b9ad928SBarry Smith EXTERN_C_END 1434b9ad928SBarry Smith 1444b9ad928SBarry Smith EXTERN_C_BEGIN 1454b9ad928SBarry Smith #undef __FUNCT__ 1464b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR" 1474b9ad928SBarry Smith int PCSORSetIterations_SOR(PC pc,int its,int lits) 1484b9ad928SBarry Smith { 1494b9ad928SBarry Smith PC_SOR *jac; 1504b9ad928SBarry Smith 1514b9ad928SBarry Smith PetscFunctionBegin; 1524b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 1534b9ad928SBarry Smith jac->its = its; 1544b9ad928SBarry Smith jac->lits = lits; 1554b9ad928SBarry Smith PetscFunctionReturn(0); 1564b9ad928SBarry Smith } 1574b9ad928SBarry Smith EXTERN_C_END 1584b9ad928SBarry Smith 1594b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1604b9ad928SBarry Smith #undef __FUNCT__ 1614b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric" 1624b9ad928SBarry Smith /*@ 1634b9ad928SBarry Smith PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR), 1644b9ad928SBarry Smith backward, or forward relaxation. The local variants perform SOR on 1654b9ad928SBarry Smith each processor. By default forward relaxation is used. 1664b9ad928SBarry Smith 1674b9ad928SBarry Smith Collective on PC 1684b9ad928SBarry Smith 1694b9ad928SBarry Smith Input Parameters: 1704b9ad928SBarry Smith + pc - the preconditioner context 1714b9ad928SBarry Smith - flag - one of the following 1724b9ad928SBarry Smith .vb 1734b9ad928SBarry Smith SOR_FORWARD_SWEEP 1744b9ad928SBarry Smith SOR_BACKWARD_SWEEP 1754b9ad928SBarry Smith SOR_SYMMETRIC_SWEEP 1764b9ad928SBarry Smith SOR_LOCAL_FORWARD_SWEEP 1774b9ad928SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 1784b9ad928SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 1794b9ad928SBarry Smith .ve 1804b9ad928SBarry Smith 1814b9ad928SBarry Smith Options Database Keys: 1824b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 1834b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 1844b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 1854b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 1864b9ad928SBarry Smith - -pc_sor_local_backward - Activates local backward version 1874b9ad928SBarry Smith 1884b9ad928SBarry Smith Notes: 1894b9ad928SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 1904b9ad928SBarry Smith which can be chosen with the option 1914b9ad928SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 1924b9ad928SBarry Smith 1934b9ad928SBarry Smith Level: intermediate 1944b9ad928SBarry Smith 1954b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric 1964b9ad928SBarry Smith 1974b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega() 1984b9ad928SBarry Smith @*/ 1994b9ad928SBarry Smith int PCSORSetSymmetric(PC pc,MatSORType flag) 2004b9ad928SBarry Smith { 2014b9ad928SBarry Smith int ierr,(*f)(PC,MatSORType); 2024b9ad928SBarry Smith 2034b9ad928SBarry Smith PetscFunctionBegin; 2044b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 2054b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetSymmetric_C",(void (**)(void))&f);CHKERRQ(ierr); 2064b9ad928SBarry Smith if (f) { 2074b9ad928SBarry Smith ierr = (*f)(pc,flag);CHKERRQ(ierr); 2084b9ad928SBarry Smith } 2094b9ad928SBarry Smith PetscFunctionReturn(0); 2104b9ad928SBarry Smith } 2114b9ad928SBarry Smith 2124b9ad928SBarry Smith #undef __FUNCT__ 2134b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega" 2144b9ad928SBarry Smith /*@ 2154b9ad928SBarry Smith PCSORSetOmega - Sets the SOR relaxation coefficient, omega 2164b9ad928SBarry Smith (where omega = 1.0 by default). 2174b9ad928SBarry Smith 2184b9ad928SBarry Smith Collective on PC 2194b9ad928SBarry Smith 2204b9ad928SBarry Smith Input Parameters: 2214b9ad928SBarry Smith + pc - the preconditioner context 2224b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2). 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith Options Database Key: 2254b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 2264b9ad928SBarry Smith 2274b9ad928SBarry Smith Level: intermediate 2284b9ad928SBarry Smith 2294b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega 2304b9ad928SBarry Smith 2314b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega() 2324b9ad928SBarry Smith @*/ 2334b9ad928SBarry Smith int PCSORSetOmega(PC pc,PetscReal omega) 2344b9ad928SBarry Smith { 2354b9ad928SBarry Smith int ierr,(*f)(PC,PetscReal); 2364b9ad928SBarry Smith 2374b9ad928SBarry Smith PetscFunctionBegin; 2384b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr); 2394b9ad928SBarry Smith if (f) { 2404b9ad928SBarry Smith ierr = (*f)(pc,omega);CHKERRQ(ierr); 2414b9ad928SBarry Smith } 2424b9ad928SBarry Smith PetscFunctionReturn(0); 2434b9ad928SBarry Smith } 2444b9ad928SBarry Smith 2454b9ad928SBarry Smith #undef __FUNCT__ 2464b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations" 2474b9ad928SBarry Smith /*@ 2484b9ad928SBarry Smith PCSORSetIterations - Sets the number of inner iterations to 2494b9ad928SBarry Smith be used by the SOR preconditioner. The default is 1. 2504b9ad928SBarry Smith 2514b9ad928SBarry Smith Collective on PC 2524b9ad928SBarry Smith 2534b9ad928SBarry Smith Input Parameters: 2544b9ad928SBarry Smith + pc - the preconditioner context 2554b9ad928SBarry Smith . lits - number of local iterations, smoothings over just variables on processor 2564b9ad928SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 2574b9ad928SBarry Smith 2584b9ad928SBarry Smith Options Database Key: 2594b9ad928SBarry Smith + -pc_sor_its <its> - Sets number of iterations 2604b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 2614b9ad928SBarry Smith 2624b9ad928SBarry Smith Level: intermediate 2634b9ad928SBarry Smith 2644b9ad928SBarry Smith Notes: When run on one processor the number of smoothings is lits*its 2654b9ad928SBarry Smith 2664b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric() 2694b9ad928SBarry Smith @*/ 2704b9ad928SBarry Smith int PCSORSetIterations(PC pc,int its,int lits) 2714b9ad928SBarry Smith { 2724b9ad928SBarry Smith int ierr,(*f)(PC,int,int); 2734b9ad928SBarry Smith 2744b9ad928SBarry Smith PetscFunctionBegin; 2754b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 2764b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetIterations_C",(void (**)(void))&f);CHKERRQ(ierr); 2774b9ad928SBarry Smith if (f) { 2784b9ad928SBarry Smith ierr = (*f)(pc,its,lits);CHKERRQ(ierr); 2794b9ad928SBarry Smith } 2804b9ad928SBarry Smith PetscFunctionReturn(0); 2814b9ad928SBarry Smith } 2824b9ad928SBarry Smith 2834b9ad928SBarry Smith /*MC 2844b9ad928SBarry Smith PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning 2854b9ad928SBarry Smith 2864b9ad928SBarry Smith Options Database Keys: 2874b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 2884b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 2894b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 2904b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 2914b9ad928SBarry Smith . -pc_sor_local_backward - Activates local backward version 2924b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 2934b9ad928SBarry Smith . -pc_sor_its <its> - Sets number of iterations 2944b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 2954b9ad928SBarry Smith 2964b9ad928SBarry Smith Level: beginner 2974b9ad928SBarry Smith 2984b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel 2994b9ad928SBarry Smith 3004b9ad928SBarry Smith Notes: Only implemented for the AIJ matrix format. 3014b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 3024b9ad928SBarry Smith Jacobi with SOR on each block. 3034b9ad928SBarry Smith 3044b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 3054b9ad928SBarry Smith PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT 3064b9ad928SBarry Smith M*/ 3074b9ad928SBarry Smith 3084b9ad928SBarry Smith EXTERN_C_BEGIN 3094b9ad928SBarry Smith #undef __FUNCT__ 3104b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR" 3114b9ad928SBarry Smith int PCCreate_SOR(PC pc) 3124b9ad928SBarry Smith { 3134b9ad928SBarry Smith int ierr; 3144b9ad928SBarry Smith PC_SOR *jac; 3154b9ad928SBarry Smith 3164b9ad928SBarry Smith PetscFunctionBegin; 3174b9ad928SBarry Smith ierr = PetscNew(PC_SOR,&jac);CHKERRQ(ierr); 3184b9ad928SBarry Smith PetscLogObjectMemory(pc,sizeof(PC_SOR)); 3194b9ad928SBarry Smith 3204b9ad928SBarry Smith pc->ops->apply = PCApply_SOR; 3214b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_SOR; 3224b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SOR; 3234b9ad928SBarry Smith pc->ops->setup = 0; 3244b9ad928SBarry Smith pc->ops->view = PCView_SOR; 3254b9ad928SBarry Smith pc->ops->destroy = PCDestroy_SOR; 3264b9ad928SBarry Smith pc->data = (void*)jac; 3274b9ad928SBarry Smith jac->sym = SOR_FORWARD_SWEEP; 3284b9ad928SBarry Smith jac->omega = 1.0; 3294b9ad928SBarry Smith jac->its = 1; 3304b9ad928SBarry Smith jac->lits = 1; 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR", 3334b9ad928SBarry Smith PCSORSetSymmetric_SOR);CHKERRQ(ierr); 3344b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR", 3354b9ad928SBarry Smith PCSORSetOmega_SOR);CHKERRQ(ierr); 3364b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR", 3374b9ad928SBarry Smith PCSORSetIterations_SOR);CHKERRQ(ierr); 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith PetscFunctionReturn(0); 3404b9ad928SBarry Smith } 3414b9ad928SBarry Smith EXTERN_C_END 3424b9ad928SBarry Smith 3434b9ad928SBarry Smith 3444b9ad928SBarry Smith 3454b9ad928SBarry Smith 3464b9ad928SBarry Smith 347