1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines a (S)SOR preconditioner for any Mat implementation 44b9ad928SBarry Smith */ 5af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 64b9ad928SBarry Smith 74b9ad928SBarry Smith typedef struct { 8c1ac3661SBarry Smith PetscInt its; /* inner iterations, number of sweeps */ 9c1ac3661SBarry Smith PetscInt lits; /* local inner iterations, number of sweeps applied by the local matrix mat->A */ 104b9ad928SBarry Smith MatSORType sym; /* forward, reverse, symmetric etc. */ 114b9ad928SBarry Smith PetscReal omega; 1229c1d7e0SHong Zhang PetscReal fshift; 134b9ad928SBarry Smith } PC_SOR; 144b9ad928SBarry Smith 154b9ad928SBarry Smith #undef __FUNCT__ 164b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_SOR" 176849ba73SBarry Smith static PetscErrorCode PCDestroy_SOR(PC pc) 184b9ad928SBarry Smith { 19dfbe8321SBarry Smith PetscErrorCode ierr; 204b9ad928SBarry Smith 214b9ad928SBarry Smith PetscFunctionBegin; 22c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 234b9ad928SBarry Smith PetscFunctionReturn(0); 244b9ad928SBarry Smith } 254b9ad928SBarry Smith 264b9ad928SBarry Smith #undef __FUNCT__ 274b9ad928SBarry Smith #define __FUNCT__ "PCApply_SOR" 286849ba73SBarry Smith static PetscErrorCode PCApply_SOR(PC pc,Vec x,Vec y) 294b9ad928SBarry Smith { 304b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 31dfbe8321SBarry Smith PetscErrorCode ierr; 32c1ac3661SBarry Smith PetscInt flag = jac->sym | SOR_ZERO_INITIAL_GUESS; 333060034bSBarry Smith PetscReal fshift; 344b9ad928SBarry Smith 354b9ad928SBarry Smith PetscFunctionBegin; 36*b0e188deSSatish Balay fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0); 37422a814eSBarry Smith ierr = MatSOR(pc->pmat,x,jac->omega,(MatSORType)flag,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" 43ace3abfcSBarry 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) 444b9ad928SBarry Smith { 454b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 46dfbe8321SBarry Smith PetscErrorCode ierr; 477319c654SBarry Smith MatSORType stype = jac->sym; 483060034bSBarry Smith PetscReal fshift; 494b9ad928SBarry Smith 504b9ad928SBarry Smith PetscFunctionBegin; 51ae15b995SBarry Smith ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr); 522fa5cd67SKarl Rupp if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS); 53*b0e188deSSatish Balay fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0); 54422a814eSBarry Smith ierr = MatSOR(pc->pmat,b,jac->omega,stype,fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr); 554d0a8057SBarry Smith *outits = its; 564d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ITS; 574b9ad928SBarry Smith PetscFunctionReturn(0); 584b9ad928SBarry Smith } 594b9ad928SBarry Smith 604b9ad928SBarry Smith #undef __FUNCT__ 614b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR" 628c34d3f5SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PetscOptions *PetscOptionsObject,PC pc) 634b9ad928SBarry Smith { 644b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 65dfbe8321SBarry Smith PetscErrorCode ierr; 66ace3abfcSBarry Smith PetscBool flg; 674b9ad928SBarry Smith 684b9ad928SBarry Smith PetscFunctionBegin; 69e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"(S)SOR options");CHKERRQ(ierr); 7094ae4db5SBarry Smith ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr); 7194ae4db5SBarry Smith ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr); 7294ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr); 7394ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr); 74acfcf0e5SJed Brown ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 754b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 76acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 774b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);} 78acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 79a9510f2eSBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);} 80acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 814b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 82acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 834b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);} 84acfcf0e5SJed Brown ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 854b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);} 864b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 874b9ad928SBarry Smith PetscFunctionReturn(0); 884b9ad928SBarry Smith } 894b9ad928SBarry Smith 904b9ad928SBarry Smith #undef __FUNCT__ 914b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR" 92dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer) 934b9ad928SBarry Smith { 944b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 954b9ad928SBarry Smith MatSORType sym = jac->sym; 962fc52814SBarry Smith const char *sortype; 97dfbe8321SBarry Smith PetscErrorCode ierr; 98ace3abfcSBarry Smith PetscBool iascii; 994b9ad928SBarry Smith 1004b9ad928SBarry Smith PetscFunctionBegin; 101251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 10232077d6dSBarry Smith if (iascii) { 1034b9ad928SBarry Smith if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer," SOR: zero initial guess\n");CHKERRQ(ierr);} 1044b9ad928SBarry Smith if (sym == SOR_APPLY_UPPER) sortype = "apply_upper"; 1054b9ad928SBarry Smith else if (sym == SOR_APPLY_LOWER) sortype = "apply_lower"; 1064b9ad928SBarry Smith else if (sym & SOR_EISENSTAT) sortype = "Eisenstat"; 107db4deed7SKarl Rupp else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP) sortype = "symmetric"; 1084b9ad928SBarry Smith else if (sym & SOR_BACKWARD_SWEEP) sortype = "backward"; 1094b9ad928SBarry Smith else if (sym & SOR_FORWARD_SWEEP) sortype = "forward"; 110db4deed7SKarl Rupp else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric"; 1114b9ad928SBarry Smith else if (sym & SOR_LOCAL_FORWARD_SWEEP) sortype = "local_forward"; 1124b9ad928SBarry Smith else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward"; 1134b9ad928SBarry Smith else sortype = "unknown"; 11457622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SOR: type = %s, iterations = %D, local iterations = %D, omega = %g\n",sortype,jac->its,jac->lits,(double)jac->omega);CHKERRQ(ierr); 1154b9ad928SBarry Smith } 1164b9ad928SBarry Smith PetscFunctionReturn(0); 1174b9ad928SBarry Smith } 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith 1204b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1214b9ad928SBarry Smith #undef __FUNCT__ 1224b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR" 123f7a08781SBarry Smith static PetscErrorCode PCSORSetSymmetric_SOR(PC pc,MatSORType flag) 1244b9ad928SBarry Smith { 125c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 1264b9ad928SBarry Smith 1274b9ad928SBarry Smith PetscFunctionBegin; 1284b9ad928SBarry Smith jac->sym = flag; 1294b9ad928SBarry Smith PetscFunctionReturn(0); 1304b9ad928SBarry Smith } 1314b9ad928SBarry Smith 1324b9ad928SBarry Smith #undef __FUNCT__ 1334b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR" 134f7a08781SBarry Smith static PetscErrorCode PCSORSetOmega_SOR(PC pc,PetscReal omega) 1354b9ad928SBarry Smith { 136c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 1374b9ad928SBarry Smith 1384b9ad928SBarry Smith PetscFunctionBegin; 139ce94432eSBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 1404b9ad928SBarry Smith jac->omega = omega; 1414b9ad928SBarry Smith PetscFunctionReturn(0); 1424b9ad928SBarry Smith } 1434b9ad928SBarry Smith 1444b9ad928SBarry Smith #undef __FUNCT__ 1454b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR" 146f7a08781SBarry Smith static PetscErrorCode PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits) 1474b9ad928SBarry Smith { 148c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 1494b9ad928SBarry Smith 1504b9ad928SBarry Smith PetscFunctionBegin; 1514b9ad928SBarry Smith jac->its = its; 1524b9ad928SBarry Smith jac->lits = lits; 1534b9ad928SBarry Smith PetscFunctionReturn(0); 1544b9ad928SBarry Smith } 1554b9ad928SBarry Smith 156c60c7ad4SBarry Smith #undef __FUNCT__ 157c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetSymmetric_SOR" 158c60c7ad4SBarry Smith static PetscErrorCode PCSORGetSymmetric_SOR(PC pc,MatSORType *flag) 159c60c7ad4SBarry Smith { 160c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 161c60c7ad4SBarry Smith 162c60c7ad4SBarry Smith PetscFunctionBegin; 163c60c7ad4SBarry Smith *flag = jac->sym; 164c60c7ad4SBarry Smith PetscFunctionReturn(0); 165c60c7ad4SBarry Smith } 166c60c7ad4SBarry Smith 167c60c7ad4SBarry Smith #undef __FUNCT__ 168c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetOmega_SOR" 169c60c7ad4SBarry Smith static PetscErrorCode PCSORGetOmega_SOR(PC pc,PetscReal *omega) 170c60c7ad4SBarry Smith { 171c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 172c60c7ad4SBarry Smith 173c60c7ad4SBarry Smith PetscFunctionBegin; 174c60c7ad4SBarry Smith *omega = jac->omega; 175c60c7ad4SBarry Smith PetscFunctionReturn(0); 176c60c7ad4SBarry Smith } 177c60c7ad4SBarry Smith 178c60c7ad4SBarry Smith #undef __FUNCT__ 179c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetIterations_SOR" 180c60c7ad4SBarry Smith static PetscErrorCode PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits) 181c60c7ad4SBarry Smith { 182c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 183c60c7ad4SBarry Smith 184c60c7ad4SBarry Smith PetscFunctionBegin; 185c60c7ad4SBarry Smith if (its) *its = jac->its; 186c60c7ad4SBarry Smith if (lits) *lits = jac->lits; 187c60c7ad4SBarry Smith PetscFunctionReturn(0); 188c60c7ad4SBarry Smith } 189c60c7ad4SBarry Smith 1904b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1914b9ad928SBarry Smith #undef __FUNCT__ 192c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetSymmetric" 193c60c7ad4SBarry Smith /*@ 1941ff2113eSBarry Smith PCSORGetSymmetric - Gets the form the SOR preconditioner is using; backward, or forward relaxation. The local variants perform SOR on 195c60c7ad4SBarry Smith each processor. By default forward relaxation is used. 196c60c7ad4SBarry Smith 197c60c7ad4SBarry Smith Logically Collective on PC 198c60c7ad4SBarry Smith 199c60c7ad4SBarry Smith Input Parameter: 200c60c7ad4SBarry Smith . pc - the preconditioner context 201c60c7ad4SBarry Smith 202c60c7ad4SBarry Smith Output Parameter: 203c60c7ad4SBarry Smith . flag - one of the following 204c60c7ad4SBarry Smith .vb 205c60c7ad4SBarry Smith SOR_FORWARD_SWEEP 206c60c7ad4SBarry Smith SOR_BACKWARD_SWEEP 207c60c7ad4SBarry Smith SOR_SYMMETRIC_SWEEP 208c60c7ad4SBarry Smith SOR_LOCAL_FORWARD_SWEEP 209c60c7ad4SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 210c60c7ad4SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 211c60c7ad4SBarry Smith .ve 212c60c7ad4SBarry Smith 213c60c7ad4SBarry Smith Options Database Keys: 214c60c7ad4SBarry Smith + -pc_sor_symmetric - Activates symmetric version 215c60c7ad4SBarry Smith . -pc_sor_backward - Activates backward version 216c60c7ad4SBarry Smith . -pc_sor_local_forward - Activates local forward version 217c60c7ad4SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 218c60c7ad4SBarry Smith - -pc_sor_local_backward - Activates local backward version 219c60c7ad4SBarry Smith 220c60c7ad4SBarry Smith Notes: 221c60c7ad4SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 222c60c7ad4SBarry Smith which can be chosen with the option 223c60c7ad4SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 224c60c7ad4SBarry Smith 225c60c7ad4SBarry Smith Level: intermediate 226c60c7ad4SBarry Smith 227c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric 228c60c7ad4SBarry Smith 229c60c7ad4SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric() 230c60c7ad4SBarry Smith @*/ 231c60c7ad4SBarry Smith PetscErrorCode PCSORGetSymmetric(PC pc,MatSORType *flag) 232c60c7ad4SBarry Smith { 233c60c7ad4SBarry Smith PetscErrorCode ierr; 234c60c7ad4SBarry Smith 235c60c7ad4SBarry Smith PetscFunctionBegin; 236c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 237c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr); 238c60c7ad4SBarry Smith PetscFunctionReturn(0); 239c60c7ad4SBarry Smith } 240c60c7ad4SBarry Smith 241c60c7ad4SBarry Smith #undef __FUNCT__ 242c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetOmega" 243c60c7ad4SBarry Smith /*@ 244c60c7ad4SBarry Smith PCSORGetOmega - Gets the SOR relaxation coefficient, omega 245c60c7ad4SBarry Smith (where omega = 1.0 by default). 246c60c7ad4SBarry Smith 247c60c7ad4SBarry Smith Logically Collective on PC 248c60c7ad4SBarry Smith 249c60c7ad4SBarry Smith Input Parameter: 250c60c7ad4SBarry Smith . pc - the preconditioner context 251c60c7ad4SBarry Smith 252c60c7ad4SBarry Smith Output Parameter: 253c60c7ad4SBarry Smith . omega - relaxation coefficient (0 < omega < 2). 254c60c7ad4SBarry Smith 255c60c7ad4SBarry Smith Options Database Key: 256c60c7ad4SBarry Smith . -pc_sor_omega <omega> - Sets omega 257c60c7ad4SBarry Smith 258c60c7ad4SBarry Smith Level: intermediate 259c60c7ad4SBarry Smith 260c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega 261c60c7ad4SBarry Smith 262c60c7ad4SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega() 263c60c7ad4SBarry Smith @*/ 264c60c7ad4SBarry Smith PetscErrorCode PCSORGetOmega(PC pc,PetscReal *omega) 265c60c7ad4SBarry Smith { 266c60c7ad4SBarry Smith PetscErrorCode ierr; 267c60c7ad4SBarry Smith 268c60c7ad4SBarry Smith PetscFunctionBegin; 269c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 270c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr); 271c60c7ad4SBarry Smith PetscFunctionReturn(0); 272c60c7ad4SBarry Smith } 273c60c7ad4SBarry Smith 274c60c7ad4SBarry Smith #undef __FUNCT__ 275c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetIterations" 276c60c7ad4SBarry Smith /*@ 277c60c7ad4SBarry Smith PCSORGetIterations - Gets the number of inner iterations to 278c60c7ad4SBarry Smith be used by the SOR preconditioner. The default is 1. 279c60c7ad4SBarry Smith 280c60c7ad4SBarry Smith Logically Collective on PC 281c60c7ad4SBarry Smith 282c60c7ad4SBarry Smith Input Parameter: 283c60c7ad4SBarry Smith . pc - the preconditioner context 284c60c7ad4SBarry Smith 285c60c7ad4SBarry Smith Output Parameter: 286c60c7ad4SBarry Smith + lits - number of local iterations, smoothings over just variables on processor 287c60c7ad4SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 288c60c7ad4SBarry Smith 289c60c7ad4SBarry Smith Options Database Key: 290c60c7ad4SBarry Smith + -pc_sor_its <its> - Sets number of iterations 291c60c7ad4SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 292c60c7ad4SBarry Smith 293c60c7ad4SBarry Smith Level: intermediate 294c60c7ad4SBarry Smith 295c60c7ad4SBarry Smith Notes: When run on one processor the number of smoothings is lits*its 296c60c7ad4SBarry Smith 297c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, iterations 298c60c7ad4SBarry Smith 299c60c7ad4SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations() 300c60c7ad4SBarry Smith @*/ 301c60c7ad4SBarry Smith PetscErrorCode PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits) 302c60c7ad4SBarry Smith { 303c60c7ad4SBarry Smith PetscErrorCode ierr; 304c60c7ad4SBarry Smith 305c60c7ad4SBarry Smith PetscFunctionBegin; 306c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 307c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr); 308c60c7ad4SBarry Smith PetscFunctionReturn(0); 309c60c7ad4SBarry Smith } 310c60c7ad4SBarry Smith 311c60c7ad4SBarry Smith #undef __FUNCT__ 3124b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric" 3134b9ad928SBarry Smith /*@ 3144b9ad928SBarry Smith PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR), 3154b9ad928SBarry Smith backward, or forward relaxation. The local variants perform SOR on 3164b9ad928SBarry Smith each processor. By default forward relaxation is used. 3174b9ad928SBarry Smith 3183f9fe445SBarry Smith Logically Collective on PC 3194b9ad928SBarry Smith 3204b9ad928SBarry Smith Input Parameters: 3214b9ad928SBarry Smith + pc - the preconditioner context 3224b9ad928SBarry Smith - flag - one of the following 3234b9ad928SBarry Smith .vb 3244b9ad928SBarry Smith SOR_FORWARD_SWEEP 3254b9ad928SBarry Smith SOR_BACKWARD_SWEEP 3264b9ad928SBarry Smith SOR_SYMMETRIC_SWEEP 3274b9ad928SBarry Smith SOR_LOCAL_FORWARD_SWEEP 3284b9ad928SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 3294b9ad928SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 3304b9ad928SBarry Smith .ve 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith Options Database Keys: 3334b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 3344b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 3354b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 3364b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 3374b9ad928SBarry Smith - -pc_sor_local_backward - Activates local backward version 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith Notes: 3404b9ad928SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 3414b9ad928SBarry Smith which can be chosen with the option 3424b9ad928SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 3434b9ad928SBarry Smith 3444b9ad928SBarry Smith Level: intermediate 3454b9ad928SBarry Smith 3464b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric 3474b9ad928SBarry Smith 3484b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega() 3494b9ad928SBarry Smith @*/ 3507087cfbeSBarry Smith PetscErrorCode PCSORSetSymmetric(PC pc,MatSORType flag) 3514b9ad928SBarry Smith { 3524ac538c5SBarry Smith PetscErrorCode ierr; 3534b9ad928SBarry Smith 3544b9ad928SBarry Smith PetscFunctionBegin; 3550700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 356c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,flag,2); 3574ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr); 3584b9ad928SBarry Smith PetscFunctionReturn(0); 3594b9ad928SBarry Smith } 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith #undef __FUNCT__ 3624b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega" 3634b9ad928SBarry Smith /*@ 3644b9ad928SBarry Smith PCSORSetOmega - Sets the SOR relaxation coefficient, omega 3654b9ad928SBarry Smith (where omega = 1.0 by default). 3664b9ad928SBarry Smith 3673f9fe445SBarry Smith Logically Collective on PC 3684b9ad928SBarry Smith 3694b9ad928SBarry Smith Input Parameters: 3704b9ad928SBarry Smith + pc - the preconditioner context 3714b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2). 3724b9ad928SBarry Smith 3734b9ad928SBarry Smith Options Database Key: 3744b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 3754b9ad928SBarry Smith 3764b9ad928SBarry Smith Level: intermediate 3774b9ad928SBarry Smith 3784b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega 3794b9ad928SBarry Smith 3804b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega() 3814b9ad928SBarry Smith @*/ 3827087cfbeSBarry Smith PetscErrorCode PCSORSetOmega(PC pc,PetscReal omega) 3834b9ad928SBarry Smith { 3844ac538c5SBarry Smith PetscErrorCode ierr; 3854b9ad928SBarry Smith 3864b9ad928SBarry Smith PetscFunctionBegin; 387c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 388c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc,omega,2); 3894ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr); 3904b9ad928SBarry Smith PetscFunctionReturn(0); 3914b9ad928SBarry Smith } 3924b9ad928SBarry Smith 3934b9ad928SBarry Smith #undef __FUNCT__ 3944b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations" 3954b9ad928SBarry Smith /*@ 3964b9ad928SBarry Smith PCSORSetIterations - Sets the number of inner iterations to 3974b9ad928SBarry Smith be used by the SOR preconditioner. The default is 1. 3984b9ad928SBarry Smith 3993f9fe445SBarry Smith Logically Collective on PC 4004b9ad928SBarry Smith 4014b9ad928SBarry Smith Input Parameters: 4024b9ad928SBarry Smith + pc - the preconditioner context 4034b9ad928SBarry Smith . lits - number of local iterations, smoothings over just variables on processor 4044b9ad928SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 4054b9ad928SBarry Smith 4064b9ad928SBarry Smith Options Database Key: 4074b9ad928SBarry Smith + -pc_sor_its <its> - Sets number of iterations 4084b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 4094b9ad928SBarry Smith 4104b9ad928SBarry Smith Level: intermediate 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith Notes: When run on one processor the number of smoothings is lits*its 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations 4154b9ad928SBarry Smith 4164b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric() 4174b9ad928SBarry Smith @*/ 4187087cfbeSBarry Smith PetscErrorCode PCSORSetIterations(PC pc,PetscInt its,PetscInt lits) 4194b9ad928SBarry Smith { 4204ac538c5SBarry Smith PetscErrorCode ierr; 4214b9ad928SBarry Smith 4224b9ad928SBarry Smith PetscFunctionBegin; 4230700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 424c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,its,2); 4254ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr); 4264b9ad928SBarry Smith PetscFunctionReturn(0); 4274b9ad928SBarry Smith } 4284b9ad928SBarry Smith 4294b9ad928SBarry Smith /*MC 4304b9ad928SBarry Smith PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning 4314b9ad928SBarry Smith 4324b9ad928SBarry Smith Options Database Keys: 4334b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 4344b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 435a9510f2eSBarry Smith . -pc_sor_forward - Activates forward version 4364b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 437a9510f2eSBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version (default version) 4384b9ad928SBarry Smith . -pc_sor_local_backward - Activates local backward version 4394b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 440422a814eSBarry Smith . -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal 441a9510f2eSBarry Smith . -pc_sor_its <its> - Sets number of iterations (default 1) 442a9510f2eSBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations (default 1) 4434b9ad928SBarry Smith 4444b9ad928SBarry Smith Level: beginner 4454b9ad928SBarry Smith 4464b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel 4474b9ad928SBarry Smith 44837a17b4dSBarry Smith Notes: Only implemented for the AIJ and SeqBAIJ matrix formats. 4494b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 4504b9ad928SBarry Smith Jacobi with SOR on each block. 4514b9ad928SBarry Smith 452422a814eSBarry Smith For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that 453422a814eSBarry Smith zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option 454422a814eSBarry Smith KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the 455422a814eSBarry Smith zero pivot. 456422a814eSBarry Smith 45737a17b4dSBarry Smith For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported. 45837a17b4dSBarry Smith 459422a814eSBarry Smith For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected 460422a814eSBarry Smith the computation is stopped with an error 461422a814eSBarry Smith 462422a814eSBarry Smith Developer Notes: We should add support for diagonal blocks that are singular to generate a Inf and thus cause KSPSolve() 463422a814eSBarry Smith to terminate with KSP_DIVERGED_NANORIF instead of stopping the program allowing 464422a814eSBarry Smith a nonlinear solver/ODE integrator to recover without stopping the program as currently happens. 465422a814eSBarry Smith 4664b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 4674b9ad928SBarry Smith PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT 4684b9ad928SBarry Smith M*/ 4694b9ad928SBarry Smith 4704b9ad928SBarry Smith #undef __FUNCT__ 4714b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR" 4728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc) 4734b9ad928SBarry Smith { 474dfbe8321SBarry Smith PetscErrorCode ierr; 4754b9ad928SBarry Smith PC_SOR *jac; 4764b9ad928SBarry Smith 4774b9ad928SBarry Smith PetscFunctionBegin; 478b00a9115SJed Brown ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 4794b9ad928SBarry Smith 4804b9ad928SBarry Smith pc->ops->apply = PCApply_SOR; 4814b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_SOR; 4824b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SOR; 4834b9ad928SBarry Smith pc->ops->setup = 0; 4844b9ad928SBarry Smith pc->ops->view = PCView_SOR; 4854b9ad928SBarry Smith pc->ops->destroy = PCDestroy_SOR; 4864b9ad928SBarry Smith pc->data = (void*)jac; 487d9bc8e36SBarry Smith jac->sym = SOR_LOCAL_SYMMETRIC_SWEEP; 4884b9ad928SBarry Smith jac->omega = 1.0; 48996fc60bcSBarry Smith jac->fshift = 0.0; 4904b9ad928SBarry Smith jac->its = 1; 4914b9ad928SBarry Smith jac->lits = 1; 4924b9ad928SBarry Smith 493bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr); 494bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr); 495bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr); 496c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr); 497c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr); 498c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr); 4994b9ad928SBarry Smith PetscFunctionReturn(0); 5004b9ad928SBarry Smith } 5014b9ad928SBarry Smith 5024b9ad928SBarry Smith 5034b9ad928SBarry Smith 5044b9ad928SBarry Smith 5054b9ad928SBarry Smith 506