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