1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 34b9ad928SBarry Smith /* 44b9ad928SBarry Smith Defines a (S)SOR preconditioner for any Mat implementation 54b9ad928SBarry Smith */ 66356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 74b9ad928SBarry Smith 84b9ad928SBarry Smith typedef struct { 9c1ac3661SBarry Smith PetscInt its; /* inner iterations, number of sweeps */ 10c1ac3661SBarry Smith PetscInt 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; 1329c1d7e0SHong Zhang PetscReal fshift; 144b9ad928SBarry Smith } PC_SOR; 154b9ad928SBarry Smith 164b9ad928SBarry Smith #undef __FUNCT__ 174b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_SOR" 186849ba73SBarry Smith static PetscErrorCode PCDestroy_SOR(PC pc) 194b9ad928SBarry Smith { 204b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 21dfbe8321SBarry Smith PetscErrorCode ierr; 224b9ad928SBarry Smith 234b9ad928SBarry Smith PetscFunctionBegin; 244b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 254b9ad928SBarry Smith PetscFunctionReturn(0); 264b9ad928SBarry Smith } 274b9ad928SBarry Smith 284b9ad928SBarry Smith #undef __FUNCT__ 294b9ad928SBarry Smith #define __FUNCT__ "PCApply_SOR" 306849ba73SBarry Smith static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y) 314b9ad928SBarry Smith { 324b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 33dfbe8321SBarry Smith PetscErrorCode ierr; 34c1ac3661SBarry Smith PetscInt flag = jac->sym | SOR_ZERO_INITIAL_GUESS; 354b9ad928SBarry Smith 364b9ad928SBarry Smith PetscFunctionBegin; 3796fc60bcSBarry Smith ierr = MatRelax(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr); 384b9ad928SBarry Smith PetscFunctionReturn(0); 394b9ad928SBarry Smith } 404b9ad928SBarry Smith 414b9ad928SBarry Smith #undef __FUNCT__ 424b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_SOR" 43*4d0a8057SBarry Smith static PetscErrorCode PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscInt *outits,PCRichardsonConvergedReason *reason) 444b9ad928SBarry Smith { 454b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 46dfbe8321SBarry Smith PetscErrorCode ierr; 474b9ad928SBarry Smith 484b9ad928SBarry Smith PetscFunctionBegin; 49ae15b995SBarry Smith ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr); 50*4d0a8057SBarry Smith ierr = MatRelax(pc->pmat,b,jac->omega,(MatSORType)jac->sym,jac->fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr); 51*4d0a8057SBarry Smith *outits = its; 52*4d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ITS; 534b9ad928SBarry Smith PetscFunctionReturn(0); 544b9ad928SBarry Smith } 554b9ad928SBarry Smith 564b9ad928SBarry Smith #undef __FUNCT__ 574b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR" 58dfbe8321SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PC pc) 594b9ad928SBarry Smith { 604b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 61dfbe8321SBarry Smith PetscErrorCode ierr; 624b9ad928SBarry Smith PetscTruth flg; 634b9ad928SBarry Smith 644b9ad928SBarry Smith PetscFunctionBegin; 654b9ad928SBarry Smith ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr); 664b9ad928SBarry Smith ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);CHKERRQ(ierr); 6796fc60bcSBarry Smith ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,0);CHKERRQ(ierr); 684b9ad928SBarry Smith ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr); 694b9ad928SBarry Smith ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);CHKERRQ(ierr); 704dc4c822SBarry Smith ierr = PetscOptionsTruthGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 714b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 724dc4c822SBarry Smith ierr = PetscOptionsTruthGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 734b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);} 74a9510f2eSBarry Smith ierr = PetscOptionsTruthGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 75a9510f2eSBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);} 76a8c7a070SBarry Smith ierr = PetscOptionsTruthGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 774b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 784dc4c822SBarry Smith ierr = PetscOptionsTruthGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 794b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);} 804dc4c822SBarry Smith ierr = PetscOptionsTruthGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 814b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);} 824b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 834b9ad928SBarry Smith PetscFunctionReturn(0); 844b9ad928SBarry Smith } 854b9ad928SBarry Smith 864b9ad928SBarry Smith #undef __FUNCT__ 874b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR" 88dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer) 894b9ad928SBarry Smith { 904b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 914b9ad928SBarry Smith MatSORType sym = jac->sym; 922fc52814SBarry Smith const char *sortype; 93dfbe8321SBarry Smith PetscErrorCode ierr; 9432077d6dSBarry Smith PetscTruth iascii; 954b9ad928SBarry Smith 964b9ad928SBarry Smith PetscFunctionBegin; 9732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 9832077d6dSBarry Smith if (iascii) { 994b9ad928SBarry Smith if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer," SOR: zero initial guess\n");CHKERRQ(ierr);} 1004b9ad928SBarry Smith if (sym == SOR_APPLY_UPPER) sortype = "apply_upper"; 1014b9ad928SBarry Smith else if (sym == SOR_APPLY_LOWER) sortype = "apply_lower"; 1024b9ad928SBarry Smith else if (sym & SOR_EISENSTAT) sortype = "Eisenstat"; 1034b9ad928SBarry Smith else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP) 1044b9ad928SBarry Smith sortype = "symmetric"; 1054b9ad928SBarry Smith else if (sym & SOR_BACKWARD_SWEEP) sortype = "backward"; 1064b9ad928SBarry Smith else if (sym & SOR_FORWARD_SWEEP) sortype = "forward"; 1074b9ad928SBarry Smith else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) 1084b9ad928SBarry Smith sortype = "local_symmetric"; 1094b9ad928SBarry Smith else if (sym & SOR_LOCAL_FORWARD_SWEEP) sortype = "local_forward"; 1104b9ad928SBarry Smith else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward"; 1114b9ad928SBarry Smith else sortype = "unknown"; 112cf576665SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SOR: type = %s, iterations = %D, local iterations = %D, omega = %G\n",sortype,jac->its,jac->lits,jac->omega);CHKERRQ(ierr); 1134b9ad928SBarry Smith } else { 11479a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name); 1154b9ad928SBarry Smith } 1164b9ad928SBarry Smith PetscFunctionReturn(0); 1174b9ad928SBarry Smith } 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith 1204b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1214b9ad928SBarry Smith EXTERN_C_BEGIN 1224b9ad928SBarry Smith #undef __FUNCT__ 1234b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR" 124dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric_SOR(PC pc,MatSORType flag) 1254b9ad928SBarry Smith { 1264b9ad928SBarry Smith PC_SOR *jac; 1274b9ad928SBarry Smith 1284b9ad928SBarry Smith PetscFunctionBegin; 1294b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 1304b9ad928SBarry Smith jac->sym = flag; 1314b9ad928SBarry Smith PetscFunctionReturn(0); 1324b9ad928SBarry Smith } 1334b9ad928SBarry Smith EXTERN_C_END 1344b9ad928SBarry Smith 1354b9ad928SBarry Smith EXTERN_C_BEGIN 1364b9ad928SBarry Smith #undef __FUNCT__ 1374b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR" 138dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega_SOR(PC pc,PetscReal omega) 1394b9ad928SBarry Smith { 1404b9ad928SBarry Smith PC_SOR *jac; 1414b9ad928SBarry Smith 1424b9ad928SBarry Smith PetscFunctionBegin; 1434b9ad928SBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 1444b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 1454b9ad928SBarry Smith jac->omega = omega; 1464b9ad928SBarry Smith PetscFunctionReturn(0); 1474b9ad928SBarry Smith } 1484b9ad928SBarry Smith EXTERN_C_END 1494b9ad928SBarry Smith 1504b9ad928SBarry Smith EXTERN_C_BEGIN 1514b9ad928SBarry Smith #undef __FUNCT__ 1524b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR" 153dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits) 1544b9ad928SBarry Smith { 1554b9ad928SBarry Smith PC_SOR *jac; 1564b9ad928SBarry Smith 1574b9ad928SBarry Smith PetscFunctionBegin; 1584b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 1594b9ad928SBarry Smith jac->its = its; 1604b9ad928SBarry Smith jac->lits = lits; 1614b9ad928SBarry Smith PetscFunctionReturn(0); 1624b9ad928SBarry Smith } 1634b9ad928SBarry Smith EXTERN_C_END 1644b9ad928SBarry Smith 1654b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1664b9ad928SBarry Smith #undef __FUNCT__ 1674b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric" 1684b9ad928SBarry Smith /*@ 1694b9ad928SBarry Smith PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR), 1704b9ad928SBarry Smith backward, or forward relaxation. The local variants perform SOR on 1714b9ad928SBarry Smith each processor. By default forward relaxation is used. 1724b9ad928SBarry Smith 1734b9ad928SBarry Smith Collective on PC 1744b9ad928SBarry Smith 1754b9ad928SBarry Smith Input Parameters: 1764b9ad928SBarry Smith + pc - the preconditioner context 1774b9ad928SBarry Smith - flag - one of the following 1784b9ad928SBarry Smith .vb 1794b9ad928SBarry Smith SOR_FORWARD_SWEEP 1804b9ad928SBarry Smith SOR_BACKWARD_SWEEP 1814b9ad928SBarry Smith SOR_SYMMETRIC_SWEEP 1824b9ad928SBarry Smith SOR_LOCAL_FORWARD_SWEEP 1834b9ad928SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 1844b9ad928SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 1854b9ad928SBarry Smith .ve 1864b9ad928SBarry Smith 1874b9ad928SBarry Smith Options Database Keys: 1884b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 1894b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 1904b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 1914b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 1924b9ad928SBarry Smith - -pc_sor_local_backward - Activates local backward version 1934b9ad928SBarry Smith 1944b9ad928SBarry Smith Notes: 1954b9ad928SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 1964b9ad928SBarry Smith which can be chosen with the option 1974b9ad928SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 1984b9ad928SBarry Smith 1994b9ad928SBarry Smith Level: intermediate 2004b9ad928SBarry Smith 2014b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric 2024b9ad928SBarry Smith 2034b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega() 2044b9ad928SBarry Smith @*/ 205dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetSymmetric(PC pc,MatSORType flag) 2064b9ad928SBarry Smith { 207dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,MatSORType); 2084b9ad928SBarry Smith 2094b9ad928SBarry Smith PetscFunctionBegin; 2104482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2114b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetSymmetric_C",(void (**)(void))&f);CHKERRQ(ierr); 2124b9ad928SBarry Smith if (f) { 2134b9ad928SBarry Smith ierr = (*f)(pc,flag);CHKERRQ(ierr); 2144b9ad928SBarry Smith } 2154b9ad928SBarry Smith PetscFunctionReturn(0); 2164b9ad928SBarry Smith } 2174b9ad928SBarry Smith 2184b9ad928SBarry Smith #undef __FUNCT__ 2194b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega" 2204b9ad928SBarry Smith /*@ 2214b9ad928SBarry Smith PCSORSetOmega - Sets the SOR relaxation coefficient, omega 2224b9ad928SBarry Smith (where omega = 1.0 by default). 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith Collective on PC 2254b9ad928SBarry Smith 2264b9ad928SBarry Smith Input Parameters: 2274b9ad928SBarry Smith + pc - the preconditioner context 2284b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2). 2294b9ad928SBarry Smith 2304b9ad928SBarry Smith Options Database Key: 2314b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 2324b9ad928SBarry Smith 2334b9ad928SBarry Smith Level: intermediate 2344b9ad928SBarry Smith 2354b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega 2364b9ad928SBarry Smith 2374b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega() 2384b9ad928SBarry Smith @*/ 239dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetOmega(PC pc,PetscReal omega) 2404b9ad928SBarry Smith { 241dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscReal); 2424b9ad928SBarry Smith 2434b9ad928SBarry Smith PetscFunctionBegin; 2444b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr); 2454b9ad928SBarry Smith if (f) { 2464b9ad928SBarry Smith ierr = (*f)(pc,omega);CHKERRQ(ierr); 2474b9ad928SBarry Smith } 2484b9ad928SBarry Smith PetscFunctionReturn(0); 2494b9ad928SBarry Smith } 2504b9ad928SBarry Smith 2514b9ad928SBarry Smith #undef __FUNCT__ 2524b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations" 2534b9ad928SBarry Smith /*@ 2544b9ad928SBarry Smith PCSORSetIterations - Sets the number of inner iterations to 2554b9ad928SBarry Smith be used by the SOR preconditioner. The default is 1. 2564b9ad928SBarry Smith 2574b9ad928SBarry Smith Collective on PC 2584b9ad928SBarry Smith 2594b9ad928SBarry Smith Input Parameters: 2604b9ad928SBarry Smith + pc - the preconditioner context 2614b9ad928SBarry Smith . lits - number of local iterations, smoothings over just variables on processor 2624b9ad928SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 2634b9ad928SBarry Smith 2644b9ad928SBarry Smith Options Database Key: 2654b9ad928SBarry Smith + -pc_sor_its <its> - Sets number of iterations 2664b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith Level: intermediate 2694b9ad928SBarry Smith 2704b9ad928SBarry Smith Notes: When run on one processor the number of smoothings is lits*its 2714b9ad928SBarry Smith 2724b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations 2734b9ad928SBarry Smith 2744b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric() 2754b9ad928SBarry Smith @*/ 276dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCSORSetIterations(PC pc,PetscInt its,PetscInt lits) 2774b9ad928SBarry Smith { 278c1ac3661SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscInt,PetscInt); 2794b9ad928SBarry Smith 2804b9ad928SBarry Smith PetscFunctionBegin; 2814482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2824b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetIterations_C",(void (**)(void))&f);CHKERRQ(ierr); 2834b9ad928SBarry Smith if (f) { 2844b9ad928SBarry Smith ierr = (*f)(pc,its,lits);CHKERRQ(ierr); 2854b9ad928SBarry Smith } 2864b9ad928SBarry Smith PetscFunctionReturn(0); 2874b9ad928SBarry Smith } 2884b9ad928SBarry Smith 2894b9ad928SBarry Smith /*MC 2904b9ad928SBarry Smith PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning 2914b9ad928SBarry Smith 2924b9ad928SBarry Smith Options Database Keys: 2934b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 2944b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 295a9510f2eSBarry Smith . -pc_sor_forward - Activates forward version 2964b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 297a9510f2eSBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version (default version) 2984b9ad928SBarry Smith . -pc_sor_local_backward - Activates local backward version 2994b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 300a9510f2eSBarry Smith . -pc_sor_its <its> - Sets number of iterations (default 1) 301a9510f2eSBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations (default 1) 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith Level: beginner 3044b9ad928SBarry Smith 3054b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel 3064b9ad928SBarry Smith 30737a17b4dSBarry Smith Notes: Only implemented for the AIJ and SeqBAIJ matrix formats. 3084b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 3094b9ad928SBarry Smith Jacobi with SOR on each block. 3104b9ad928SBarry Smith 31137a17b4dSBarry Smith For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported. 31237a17b4dSBarry Smith 3134b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 3144b9ad928SBarry Smith PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT 3154b9ad928SBarry Smith M*/ 3164b9ad928SBarry Smith 3174b9ad928SBarry Smith EXTERN_C_BEGIN 3184b9ad928SBarry Smith #undef __FUNCT__ 3194b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR" 320dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_SOR(PC pc) 3214b9ad928SBarry Smith { 322dfbe8321SBarry Smith PetscErrorCode ierr; 3234b9ad928SBarry Smith PC_SOR *jac; 3244b9ad928SBarry Smith 3254b9ad928SBarry Smith PetscFunctionBegin; 32638f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_SOR,&jac);CHKERRQ(ierr); 3274b9ad928SBarry Smith 3284b9ad928SBarry Smith pc->ops->apply = PCApply_SOR; 3294b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_SOR; 3304b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SOR; 3314b9ad928SBarry Smith pc->ops->setup = 0; 3324b9ad928SBarry Smith pc->ops->view = PCView_SOR; 3334b9ad928SBarry Smith pc->ops->destroy = PCDestroy_SOR; 3344b9ad928SBarry Smith pc->data = (void*)jac; 335d9bc8e36SBarry Smith jac->sym = SOR_LOCAL_SYMMETRIC_SWEEP; 3364b9ad928SBarry Smith jac->omega = 1.0; 33796fc60bcSBarry Smith jac->fshift = 0.0; 3384b9ad928SBarry Smith jac->its = 1; 3394b9ad928SBarry Smith jac->lits = 1; 3404b9ad928SBarry Smith 3414b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR", 3424b9ad928SBarry Smith PCSORSetSymmetric_SOR);CHKERRQ(ierr); 3434b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR", 3444b9ad928SBarry Smith PCSORSetOmega_SOR);CHKERRQ(ierr); 3454b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR", 3464b9ad928SBarry Smith PCSORSetIterations_SOR);CHKERRQ(ierr); 3474b9ad928SBarry Smith 3484b9ad928SBarry Smith PetscFunctionReturn(0); 3494b9ad928SBarry Smith } 3504b9ad928SBarry Smith EXTERN_C_END 3514b9ad928SBarry Smith 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith 3544b9ad928SBarry Smith 3554b9ad928SBarry Smith 356