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 { 20dfbe8321SBarry Smith PetscErrorCode ierr; 214b9ad928SBarry Smith 224b9ad928SBarry Smith PetscFunctionBegin; 23*c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 244b9ad928SBarry Smith PetscFunctionReturn(0); 254b9ad928SBarry Smith } 264b9ad928SBarry Smith 274b9ad928SBarry Smith #undef __FUNCT__ 284b9ad928SBarry Smith #define __FUNCT__ "PCApply_SOR" 296849ba73SBarry Smith static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y) 304b9ad928SBarry Smith { 314b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 32dfbe8321SBarry Smith PetscErrorCode ierr; 33c1ac3661SBarry Smith PetscInt flag = jac->sym | SOR_ZERO_INITIAL_GUESS; 344b9ad928SBarry Smith 354b9ad928SBarry Smith PetscFunctionBegin; 3641f059aeSBarry Smith ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,jac->fshift,jac->its,jac->lits,y);CHKERRQ(ierr); 374b9ad928SBarry Smith PetscFunctionReturn(0); 384b9ad928SBarry Smith } 394b9ad928SBarry Smith 404b9ad928SBarry Smith #undef __FUNCT__ 414b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_SOR" 42ace3abfcSBarry Smith static PetscErrorCode PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool guesszero,PetscInt *outits,PCRichardsonConvergedReason *reason) 434b9ad928SBarry Smith { 444b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 45dfbe8321SBarry Smith PetscErrorCode ierr; 467319c654SBarry Smith MatSORType stype = jac->sym; 474b9ad928SBarry Smith 484b9ad928SBarry Smith PetscFunctionBegin; 49ae15b995SBarry Smith ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr); 507319c654SBarry Smith if (guesszero) { 51de8939e3SBarry Smith stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS); 527319c654SBarry Smith } 5341f059aeSBarry Smith ierr = MatSOR(pc->pmat,b,jac->omega,stype,jac->fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr); 544d0a8057SBarry Smith *outits = its; 554d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ITS; 564b9ad928SBarry Smith PetscFunctionReturn(0); 574b9ad928SBarry Smith } 584b9ad928SBarry Smith 594b9ad928SBarry Smith #undef __FUNCT__ 604b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR" 61dfbe8321SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PC pc) 624b9ad928SBarry Smith { 634b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 64dfbe8321SBarry Smith PetscErrorCode ierr; 65ace3abfcSBarry Smith PetscBool flg; 664b9ad928SBarry Smith 674b9ad928SBarry Smith PetscFunctionBegin; 684b9ad928SBarry Smith ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr); 694b9ad928SBarry Smith ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);CHKERRQ(ierr); 7096fc60bcSBarry Smith ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,0);CHKERRQ(ierr); 714b9ad928SBarry Smith ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr); 724b9ad928SBarry Smith ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);CHKERRQ(ierr); 73acfcf0e5SJed Brown ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 744b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 75acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 764b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);} 77acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 78a9510f2eSBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);} 79acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 804b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 81acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 824b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);} 83acfcf0e5SJed Brown ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 844b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);} 854b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 864b9ad928SBarry Smith PetscFunctionReturn(0); 874b9ad928SBarry Smith } 884b9ad928SBarry Smith 894b9ad928SBarry Smith #undef __FUNCT__ 904b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR" 91dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer) 924b9ad928SBarry Smith { 934b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 944b9ad928SBarry Smith MatSORType sym = jac->sym; 952fc52814SBarry Smith const char *sortype; 96dfbe8321SBarry Smith PetscErrorCode ierr; 97ace3abfcSBarry Smith PetscBool iascii; 984b9ad928SBarry Smith 994b9ad928SBarry Smith PetscFunctionBegin; 1002692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 10132077d6dSBarry Smith if (iascii) { 1024b9ad928SBarry Smith if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer," SOR: zero initial guess\n");CHKERRQ(ierr);} 1034b9ad928SBarry Smith if (sym == SOR_APPLY_UPPER) sortype = "apply_upper"; 1044b9ad928SBarry Smith else if (sym == SOR_APPLY_LOWER) sortype = "apply_lower"; 1054b9ad928SBarry Smith else if (sym & SOR_EISENSTAT) sortype = "Eisenstat"; 1064b9ad928SBarry Smith else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP) 1074b9ad928SBarry Smith sortype = "symmetric"; 1084b9ad928SBarry Smith else if (sym & SOR_BACKWARD_SWEEP) sortype = "backward"; 1094b9ad928SBarry Smith else if (sym & SOR_FORWARD_SWEEP) sortype = "forward"; 1104b9ad928SBarry Smith else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) 1114b9ad928SBarry Smith sortype = "local_symmetric"; 1124b9ad928SBarry Smith else if (sym & SOR_LOCAL_FORWARD_SWEEP) sortype = "local_forward"; 1134b9ad928SBarry Smith else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward"; 1144b9ad928SBarry Smith else sortype = "unknown"; 115cf576665SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SOR: type = %s, iterations = %D, local iterations = %D, omega = %G\n",sortype,jac->its,jac->lits,jac->omega);CHKERRQ(ierr); 1164b9ad928SBarry Smith } else { 11765e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name); 1184b9ad928SBarry Smith } 1194b9ad928SBarry Smith PetscFunctionReturn(0); 1204b9ad928SBarry Smith } 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith 1234b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1244b9ad928SBarry Smith EXTERN_C_BEGIN 1254b9ad928SBarry Smith #undef __FUNCT__ 1264b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR" 1277087cfbeSBarry Smith PetscErrorCode PCSORSetSymmetric_SOR(PC pc,MatSORType flag) 1284b9ad928SBarry Smith { 1294b9ad928SBarry Smith PC_SOR *jac; 1304b9ad928SBarry Smith 1314b9ad928SBarry Smith PetscFunctionBegin; 1324b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 1334b9ad928SBarry Smith jac->sym = flag; 1344b9ad928SBarry Smith PetscFunctionReturn(0); 1354b9ad928SBarry Smith } 1364b9ad928SBarry Smith EXTERN_C_END 1374b9ad928SBarry Smith 1384b9ad928SBarry Smith EXTERN_C_BEGIN 1394b9ad928SBarry Smith #undef __FUNCT__ 1404b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR" 1417087cfbeSBarry Smith PetscErrorCode PCSORSetOmega_SOR(PC pc,PetscReal omega) 1424b9ad928SBarry Smith { 1434b9ad928SBarry Smith PC_SOR *jac; 1444b9ad928SBarry Smith 1454b9ad928SBarry Smith PetscFunctionBegin; 146e7e72b3dSBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 1474b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 1484b9ad928SBarry Smith jac->omega = omega; 1494b9ad928SBarry Smith PetscFunctionReturn(0); 1504b9ad928SBarry Smith } 1514b9ad928SBarry Smith EXTERN_C_END 1524b9ad928SBarry Smith 1534b9ad928SBarry Smith EXTERN_C_BEGIN 1544b9ad928SBarry Smith #undef __FUNCT__ 1554b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR" 1567087cfbeSBarry Smith PetscErrorCode PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits) 1574b9ad928SBarry Smith { 1584b9ad928SBarry Smith PC_SOR *jac; 1594b9ad928SBarry Smith 1604b9ad928SBarry Smith PetscFunctionBegin; 1614b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 1624b9ad928SBarry Smith jac->its = its; 1634b9ad928SBarry Smith jac->lits = lits; 1644b9ad928SBarry Smith PetscFunctionReturn(0); 1654b9ad928SBarry Smith } 1664b9ad928SBarry Smith EXTERN_C_END 1674b9ad928SBarry Smith 1684b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1694b9ad928SBarry Smith #undef __FUNCT__ 1704b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric" 1714b9ad928SBarry Smith /*@ 1724b9ad928SBarry Smith PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR), 1734b9ad928SBarry Smith backward, or forward relaxation. The local variants perform SOR on 1744b9ad928SBarry Smith each processor. By default forward relaxation is used. 1754b9ad928SBarry Smith 1763f9fe445SBarry Smith Logically Collective on PC 1774b9ad928SBarry Smith 1784b9ad928SBarry Smith Input Parameters: 1794b9ad928SBarry Smith + pc - the preconditioner context 1804b9ad928SBarry Smith - flag - one of the following 1814b9ad928SBarry Smith .vb 1824b9ad928SBarry Smith SOR_FORWARD_SWEEP 1834b9ad928SBarry Smith SOR_BACKWARD_SWEEP 1844b9ad928SBarry Smith SOR_SYMMETRIC_SWEEP 1854b9ad928SBarry Smith SOR_LOCAL_FORWARD_SWEEP 1864b9ad928SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 1874b9ad928SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 1884b9ad928SBarry Smith .ve 1894b9ad928SBarry Smith 1904b9ad928SBarry Smith Options Database Keys: 1914b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 1924b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 1934b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 1944b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 1954b9ad928SBarry Smith - -pc_sor_local_backward - Activates local backward version 1964b9ad928SBarry Smith 1974b9ad928SBarry Smith Notes: 1984b9ad928SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 1994b9ad928SBarry Smith which can be chosen with the option 2004b9ad928SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 2014b9ad928SBarry Smith 2024b9ad928SBarry Smith Level: intermediate 2034b9ad928SBarry Smith 2044b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric 2054b9ad928SBarry Smith 2064b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega() 2074b9ad928SBarry Smith @*/ 2087087cfbeSBarry Smith PetscErrorCode PCSORSetSymmetric(PC pc,MatSORType flag) 2094b9ad928SBarry Smith { 2104ac538c5SBarry Smith PetscErrorCode ierr; 2114b9ad928SBarry Smith 2124b9ad928SBarry Smith PetscFunctionBegin; 2130700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 214c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,flag,2); 2154ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr); 2164b9ad928SBarry Smith PetscFunctionReturn(0); 2174b9ad928SBarry Smith } 2184b9ad928SBarry Smith 2194b9ad928SBarry Smith #undef __FUNCT__ 2204b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega" 2214b9ad928SBarry Smith /*@ 2224b9ad928SBarry Smith PCSORSetOmega - Sets the SOR relaxation coefficient, omega 2234b9ad928SBarry Smith (where omega = 1.0 by default). 2244b9ad928SBarry Smith 2253f9fe445SBarry Smith Logically Collective on PC 2264b9ad928SBarry Smith 2274b9ad928SBarry Smith Input Parameters: 2284b9ad928SBarry Smith + pc - the preconditioner context 2294b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2). 2304b9ad928SBarry Smith 2314b9ad928SBarry Smith Options Database Key: 2324b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 2334b9ad928SBarry Smith 2344b9ad928SBarry Smith Level: intermediate 2354b9ad928SBarry Smith 2364b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega 2374b9ad928SBarry Smith 2384b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega() 2394b9ad928SBarry Smith @*/ 2407087cfbeSBarry Smith PetscErrorCode PCSORSetOmega(PC pc,PetscReal omega) 2414b9ad928SBarry Smith { 2424ac538c5SBarry Smith PetscErrorCode ierr; 2434b9ad928SBarry Smith 2444b9ad928SBarry Smith PetscFunctionBegin; 245c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 246c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc,omega,2); 2474ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr); 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 2573f9fe445SBarry Smith Logically 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 @*/ 2767087cfbeSBarry Smith PetscErrorCode PCSORSetIterations(PC pc,PetscInt its,PetscInt lits) 2774b9ad928SBarry Smith { 2784ac538c5SBarry Smith PetscErrorCode ierr; 2794b9ad928SBarry Smith 2804b9ad928SBarry Smith PetscFunctionBegin; 2810700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 282c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,its,2); 2834ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr); 2844b9ad928SBarry Smith PetscFunctionReturn(0); 2854b9ad928SBarry Smith } 2864b9ad928SBarry Smith 2874b9ad928SBarry Smith /*MC 2884b9ad928SBarry Smith PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning 2894b9ad928SBarry Smith 2904b9ad928SBarry Smith Options Database Keys: 2914b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 2924b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 293a9510f2eSBarry Smith . -pc_sor_forward - Activates forward version 2944b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 295a9510f2eSBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version (default version) 2964b9ad928SBarry Smith . -pc_sor_local_backward - Activates local backward version 2974b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 298a9510f2eSBarry Smith . -pc_sor_its <its> - Sets number of iterations (default 1) 299a9510f2eSBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations (default 1) 3004b9ad928SBarry Smith 3014b9ad928SBarry Smith Level: beginner 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel 3044b9ad928SBarry Smith 30537a17b4dSBarry Smith Notes: Only implemented for the AIJ and SeqBAIJ matrix formats. 3064b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 3074b9ad928SBarry Smith Jacobi with SOR on each block. 3084b9ad928SBarry Smith 30937a17b4dSBarry Smith For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported. 31037a17b4dSBarry Smith 3114b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 3124b9ad928SBarry Smith PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT 3134b9ad928SBarry Smith M*/ 3144b9ad928SBarry Smith 3154b9ad928SBarry Smith EXTERN_C_BEGIN 3164b9ad928SBarry Smith #undef __FUNCT__ 3174b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR" 3187087cfbeSBarry Smith PetscErrorCode PCCreate_SOR(PC pc) 3194b9ad928SBarry Smith { 320dfbe8321SBarry Smith PetscErrorCode ierr; 3214b9ad928SBarry Smith PC_SOR *jac; 3224b9ad928SBarry Smith 3234b9ad928SBarry Smith PetscFunctionBegin; 32438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_SOR,&jac);CHKERRQ(ierr); 3254b9ad928SBarry Smith 3264b9ad928SBarry Smith pc->ops->apply = PCApply_SOR; 3274b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_SOR; 3284b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SOR; 3294b9ad928SBarry Smith pc->ops->setup = 0; 3304b9ad928SBarry Smith pc->ops->view = PCView_SOR; 3314b9ad928SBarry Smith pc->ops->destroy = PCDestroy_SOR; 3324b9ad928SBarry Smith pc->data = (void*)jac; 333d9bc8e36SBarry Smith jac->sym = SOR_LOCAL_SYMMETRIC_SWEEP; 3344b9ad928SBarry Smith jac->omega = 1.0; 33596fc60bcSBarry Smith jac->fshift = 0.0; 3364b9ad928SBarry Smith jac->its = 1; 3374b9ad928SBarry Smith jac->lits = 1; 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR", 3404b9ad928SBarry Smith PCSORSetSymmetric_SOR);CHKERRQ(ierr); 3414b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR", 3424b9ad928SBarry Smith PCSORSetOmega_SOR);CHKERRQ(ierr); 3434b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR", 3444b9ad928SBarry Smith PCSORSetIterations_SOR);CHKERRQ(ierr); 3454b9ad928SBarry Smith 3464b9ad928SBarry Smith PetscFunctionReturn(0); 3474b9ad928SBarry Smith } 3484b9ad928SBarry Smith EXTERN_C_END 3494b9ad928SBarry Smith 3504b9ad928SBarry Smith 3514b9ad928SBarry Smith 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith 354