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; 36b0e188deSSatish 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); 38*539c167fSBarry Smith ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr); 394b9ad928SBarry Smith PetscFunctionReturn(0); 404b9ad928SBarry Smith } 414b9ad928SBarry Smith 424b9ad928SBarry Smith #undef __FUNCT__ 434b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_SOR" 44ace3abfcSBarry 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) 454b9ad928SBarry Smith { 464b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 47dfbe8321SBarry Smith PetscErrorCode ierr; 487319c654SBarry Smith MatSORType stype = jac->sym; 493060034bSBarry Smith PetscReal fshift; 504b9ad928SBarry Smith 514b9ad928SBarry Smith PetscFunctionBegin; 52ae15b995SBarry Smith ierr = PetscInfo1(pc,"Warning, convergence critera ignored, using %D iterations\n",its);CHKERRQ(ierr); 532fa5cd67SKarl Rupp if (guesszero) stype = (MatSORType) (stype | SOR_ZERO_INITIAL_GUESS); 54b0e188deSSatish Balay fshift = (jac->fshift ? jac->fshift : pc->erroriffailure ? 0.0 : -1.0); 55422a814eSBarry Smith ierr = MatSOR(pc->pmat,b,jac->omega,stype,fshift,its*jac->its,jac->lits,y);CHKERRQ(ierr); 56*539c167fSBarry Smith ierr = MatFactorGetError(pc->pmat,(MatFactorError*)&pc->failedreason);CHKERRQ(ierr); 574d0a8057SBarry Smith *outits = its; 584d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ITS; 594b9ad928SBarry Smith PetscFunctionReturn(0); 604b9ad928SBarry Smith } 614b9ad928SBarry Smith 624b9ad928SBarry Smith #undef __FUNCT__ 634b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR" 644416b707SBarry Smith PetscErrorCode PCSetFromOptions_SOR(PetscOptionItems *PetscOptionsObject,PC pc) 654b9ad928SBarry Smith { 664b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 67dfbe8321SBarry Smith PetscErrorCode ierr; 68ace3abfcSBarry Smith PetscBool flg; 694b9ad928SBarry Smith 704b9ad928SBarry Smith PetscFunctionBegin; 71e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"(S)SOR options");CHKERRQ(ierr); 7294ae4db5SBarry Smith ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,NULL);CHKERRQ(ierr); 7394ae4db5SBarry Smith ierr = PetscOptionsReal("-pc_sor_diagonal_shift","Add to the diagonal entries","",jac->fshift,&jac->fshift,NULL);CHKERRQ(ierr); 7494ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,NULL);CHKERRQ(ierr); 7594ae4db5SBarry Smith ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,NULL);CHKERRQ(ierr); 76acfcf0e5SJed Brown ierr = PetscOptionsBoolGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 774b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 78acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 794b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);} 80acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_forward","use forward sweep","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 81a9510f2eSBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_FORWARD_SWEEP);CHKERRQ(ierr);} 82acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_local_symmetric","use SSOR separately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 834b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 84acfcf0e5SJed Brown ierr = PetscOptionsBoolGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 854b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);} 86acfcf0e5SJed Brown ierr = PetscOptionsBoolGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 874b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);} 884b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 894b9ad928SBarry Smith PetscFunctionReturn(0); 904b9ad928SBarry Smith } 914b9ad928SBarry Smith 924b9ad928SBarry Smith #undef __FUNCT__ 934b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR" 94dfbe8321SBarry Smith PetscErrorCode PCView_SOR(PC pc,PetscViewer viewer) 954b9ad928SBarry Smith { 964b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 974b9ad928SBarry Smith MatSORType sym = jac->sym; 982fc52814SBarry Smith const char *sortype; 99dfbe8321SBarry Smith PetscErrorCode ierr; 100ace3abfcSBarry Smith PetscBool iascii; 1014b9ad928SBarry Smith 1024b9ad928SBarry Smith PetscFunctionBegin; 103251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 10432077d6dSBarry Smith if (iascii) { 1054b9ad928SBarry Smith if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer," SOR: zero initial guess\n");CHKERRQ(ierr);} 1064b9ad928SBarry Smith if (sym == SOR_APPLY_UPPER) sortype = "apply_upper"; 1074b9ad928SBarry Smith else if (sym == SOR_APPLY_LOWER) sortype = "apply_lower"; 1084b9ad928SBarry Smith else if (sym & SOR_EISENSTAT) sortype = "Eisenstat"; 109db4deed7SKarl Rupp else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP) sortype = "symmetric"; 1104b9ad928SBarry Smith else if (sym & SOR_BACKWARD_SWEEP) sortype = "backward"; 1114b9ad928SBarry Smith else if (sym & SOR_FORWARD_SWEEP) sortype = "forward"; 112db4deed7SKarl Rupp else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric"; 1134b9ad928SBarry Smith else if (sym & SOR_LOCAL_FORWARD_SWEEP) sortype = "local_forward"; 1144b9ad928SBarry Smith else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward"; 1154b9ad928SBarry Smith else sortype = "unknown"; 11657622a8eSBarry 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); 1174b9ad928SBarry Smith } 1184b9ad928SBarry Smith PetscFunctionReturn(0); 1194b9ad928SBarry Smith } 1204b9ad928SBarry Smith 1214b9ad928SBarry Smith 1224b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1234b9ad928SBarry Smith #undef __FUNCT__ 1244b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR" 125f7a08781SBarry Smith static PetscErrorCode PCSORSetSymmetric_SOR(PC pc,MatSORType flag) 1264b9ad928SBarry Smith { 127c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith PetscFunctionBegin; 1304b9ad928SBarry Smith jac->sym = flag; 1314b9ad928SBarry Smith PetscFunctionReturn(0); 1324b9ad928SBarry Smith } 1334b9ad928SBarry Smith 1344b9ad928SBarry Smith #undef __FUNCT__ 1354b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR" 136f7a08781SBarry Smith static PetscErrorCode PCSORSetOmega_SOR(PC pc,PetscReal omega) 1374b9ad928SBarry Smith { 138c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 1394b9ad928SBarry Smith 1404b9ad928SBarry Smith PetscFunctionBegin; 141ce94432eSBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 1424b9ad928SBarry Smith jac->omega = omega; 1434b9ad928SBarry Smith PetscFunctionReturn(0); 1444b9ad928SBarry Smith } 1454b9ad928SBarry Smith 1464b9ad928SBarry Smith #undef __FUNCT__ 1474b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR" 148f7a08781SBarry Smith static PetscErrorCode PCSORSetIterations_SOR(PC pc,PetscInt its,PetscInt lits) 1494b9ad928SBarry Smith { 150c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 1514b9ad928SBarry Smith 1524b9ad928SBarry Smith PetscFunctionBegin; 1534b9ad928SBarry Smith jac->its = its; 1544b9ad928SBarry Smith jac->lits = lits; 1554b9ad928SBarry Smith PetscFunctionReturn(0); 1564b9ad928SBarry Smith } 1574b9ad928SBarry Smith 158c60c7ad4SBarry Smith #undef __FUNCT__ 159c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetSymmetric_SOR" 160c60c7ad4SBarry Smith static PetscErrorCode PCSORGetSymmetric_SOR(PC pc,MatSORType *flag) 161c60c7ad4SBarry Smith { 162c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 163c60c7ad4SBarry Smith 164c60c7ad4SBarry Smith PetscFunctionBegin; 165c60c7ad4SBarry Smith *flag = jac->sym; 166c60c7ad4SBarry Smith PetscFunctionReturn(0); 167c60c7ad4SBarry Smith } 168c60c7ad4SBarry Smith 169c60c7ad4SBarry Smith #undef __FUNCT__ 170c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetOmega_SOR" 171c60c7ad4SBarry Smith static PetscErrorCode PCSORGetOmega_SOR(PC pc,PetscReal *omega) 172c60c7ad4SBarry Smith { 173c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 174c60c7ad4SBarry Smith 175c60c7ad4SBarry Smith PetscFunctionBegin; 176c60c7ad4SBarry Smith *omega = jac->omega; 177c60c7ad4SBarry Smith PetscFunctionReturn(0); 178c60c7ad4SBarry Smith } 179c60c7ad4SBarry Smith 180c60c7ad4SBarry Smith #undef __FUNCT__ 181c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetIterations_SOR" 182c60c7ad4SBarry Smith static PetscErrorCode PCSORGetIterations_SOR(PC pc,PetscInt *its,PetscInt *lits) 183c60c7ad4SBarry Smith { 184c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 185c60c7ad4SBarry Smith 186c60c7ad4SBarry Smith PetscFunctionBegin; 187c60c7ad4SBarry Smith if (its) *its = jac->its; 188c60c7ad4SBarry Smith if (lits) *lits = jac->lits; 189c60c7ad4SBarry Smith PetscFunctionReturn(0); 190c60c7ad4SBarry Smith } 191c60c7ad4SBarry Smith 1924b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 1934b9ad928SBarry Smith #undef __FUNCT__ 194c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetSymmetric" 195c60c7ad4SBarry Smith /*@ 1961ff2113eSBarry Smith PCSORGetSymmetric - Gets the form the SOR preconditioner is using; backward, or forward relaxation. The local variants perform SOR on 197c60c7ad4SBarry Smith each processor. By default forward relaxation is used. 198c60c7ad4SBarry Smith 199c60c7ad4SBarry Smith Logically Collective on PC 200c60c7ad4SBarry Smith 201c60c7ad4SBarry Smith Input Parameter: 202c60c7ad4SBarry Smith . pc - the preconditioner context 203c60c7ad4SBarry Smith 204c60c7ad4SBarry Smith Output Parameter: 205c60c7ad4SBarry Smith . flag - one of the following 206c60c7ad4SBarry Smith .vb 207c60c7ad4SBarry Smith SOR_FORWARD_SWEEP 208c60c7ad4SBarry Smith SOR_BACKWARD_SWEEP 209c60c7ad4SBarry Smith SOR_SYMMETRIC_SWEEP 210c60c7ad4SBarry Smith SOR_LOCAL_FORWARD_SWEEP 211c60c7ad4SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 212c60c7ad4SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 213c60c7ad4SBarry Smith .ve 214c60c7ad4SBarry Smith 215c60c7ad4SBarry Smith Options Database Keys: 216c60c7ad4SBarry Smith + -pc_sor_symmetric - Activates symmetric version 217c60c7ad4SBarry Smith . -pc_sor_backward - Activates backward version 218c60c7ad4SBarry Smith . -pc_sor_local_forward - Activates local forward version 219c60c7ad4SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 220c60c7ad4SBarry Smith - -pc_sor_local_backward - Activates local backward version 221c60c7ad4SBarry Smith 222c60c7ad4SBarry Smith Notes: 223c60c7ad4SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 224c60c7ad4SBarry Smith which can be chosen with the option 225c60c7ad4SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 226c60c7ad4SBarry Smith 227c60c7ad4SBarry Smith Level: intermediate 228c60c7ad4SBarry Smith 229c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric 230c60c7ad4SBarry Smith 231c60c7ad4SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega(), PCSORSetSymmetric() 232c60c7ad4SBarry Smith @*/ 233c60c7ad4SBarry Smith PetscErrorCode PCSORGetSymmetric(PC pc,MatSORType *flag) 234c60c7ad4SBarry Smith { 235c60c7ad4SBarry Smith PetscErrorCode ierr; 236c60c7ad4SBarry Smith 237c60c7ad4SBarry Smith PetscFunctionBegin; 238c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 239c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCSORGetSymmetric_C",(PC,MatSORType*),(pc,flag));CHKERRQ(ierr); 240c60c7ad4SBarry Smith PetscFunctionReturn(0); 241c60c7ad4SBarry Smith } 242c60c7ad4SBarry Smith 243c60c7ad4SBarry Smith #undef __FUNCT__ 244c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetOmega" 245c60c7ad4SBarry Smith /*@ 246c60c7ad4SBarry Smith PCSORGetOmega - Gets the SOR relaxation coefficient, omega 247c60c7ad4SBarry Smith (where omega = 1.0 by default). 248c60c7ad4SBarry Smith 249c60c7ad4SBarry Smith Logically Collective on PC 250c60c7ad4SBarry Smith 251c60c7ad4SBarry Smith Input Parameter: 252c60c7ad4SBarry Smith . pc - the preconditioner context 253c60c7ad4SBarry Smith 254c60c7ad4SBarry Smith Output Parameter: 255c60c7ad4SBarry Smith . omega - relaxation coefficient (0 < omega < 2). 256c60c7ad4SBarry Smith 257c60c7ad4SBarry Smith Options Database Key: 258c60c7ad4SBarry Smith . -pc_sor_omega <omega> - Sets omega 259c60c7ad4SBarry Smith 260c60c7ad4SBarry Smith Level: intermediate 261c60c7ad4SBarry Smith 262c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega 263c60c7ad4SBarry Smith 264c60c7ad4SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega(), PCSORSetOmega() 265c60c7ad4SBarry Smith @*/ 266c60c7ad4SBarry Smith PetscErrorCode PCSORGetOmega(PC pc,PetscReal *omega) 267c60c7ad4SBarry Smith { 268c60c7ad4SBarry Smith PetscErrorCode ierr; 269c60c7ad4SBarry Smith 270c60c7ad4SBarry Smith PetscFunctionBegin; 271c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 272c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCSORGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr); 273c60c7ad4SBarry Smith PetscFunctionReturn(0); 274c60c7ad4SBarry Smith } 275c60c7ad4SBarry Smith 276c60c7ad4SBarry Smith #undef __FUNCT__ 277c60c7ad4SBarry Smith #define __FUNCT__ "PCSORGetIterations" 278c60c7ad4SBarry Smith /*@ 279c60c7ad4SBarry Smith PCSORGetIterations - Gets the number of inner iterations to 280c60c7ad4SBarry Smith be used by the SOR preconditioner. The default is 1. 281c60c7ad4SBarry Smith 282c60c7ad4SBarry Smith Logically Collective on PC 283c60c7ad4SBarry Smith 284c60c7ad4SBarry Smith Input Parameter: 285c60c7ad4SBarry Smith . pc - the preconditioner context 286c60c7ad4SBarry Smith 287c60c7ad4SBarry Smith Output Parameter: 288c60c7ad4SBarry Smith + lits - number of local iterations, smoothings over just variables on processor 289c60c7ad4SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 290c60c7ad4SBarry Smith 291c60c7ad4SBarry Smith Options Database Key: 292c60c7ad4SBarry Smith + -pc_sor_its <its> - Sets number of iterations 293c60c7ad4SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 294c60c7ad4SBarry Smith 295c60c7ad4SBarry Smith Level: intermediate 296c60c7ad4SBarry Smith 297c60c7ad4SBarry Smith Notes: When run on one processor the number of smoothings is lits*its 298c60c7ad4SBarry Smith 299c60c7ad4SBarry Smith .keywords: PC, SOR, SSOR, set, iterations 300c60c7ad4SBarry Smith 301c60c7ad4SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric(), PCSORSetIterations() 302c60c7ad4SBarry Smith @*/ 303c60c7ad4SBarry Smith PetscErrorCode PCSORGetIterations(PC pc,PetscInt *its,PetscInt *lits) 304c60c7ad4SBarry Smith { 305c60c7ad4SBarry Smith PetscErrorCode ierr; 306c60c7ad4SBarry Smith 307c60c7ad4SBarry Smith PetscFunctionBegin; 308c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 309c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCSORGetIterations_C",(PC,PetscInt*,PetscInt*),(pc,its,lits));CHKERRQ(ierr); 310c60c7ad4SBarry Smith PetscFunctionReturn(0); 311c60c7ad4SBarry Smith } 312c60c7ad4SBarry Smith 313c60c7ad4SBarry Smith #undef __FUNCT__ 3144b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric" 3154b9ad928SBarry Smith /*@ 3164b9ad928SBarry Smith PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR), 3174b9ad928SBarry Smith backward, or forward relaxation. The local variants perform SOR on 3184b9ad928SBarry Smith each processor. By default forward relaxation is used. 3194b9ad928SBarry Smith 3203f9fe445SBarry Smith Logically Collective on PC 3214b9ad928SBarry Smith 3224b9ad928SBarry Smith Input Parameters: 3234b9ad928SBarry Smith + pc - the preconditioner context 3244b9ad928SBarry Smith - flag - one of the following 3254b9ad928SBarry Smith .vb 3264b9ad928SBarry Smith SOR_FORWARD_SWEEP 3274b9ad928SBarry Smith SOR_BACKWARD_SWEEP 3284b9ad928SBarry Smith SOR_SYMMETRIC_SWEEP 3294b9ad928SBarry Smith SOR_LOCAL_FORWARD_SWEEP 3304b9ad928SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 3314b9ad928SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 3324b9ad928SBarry Smith .ve 3334b9ad928SBarry Smith 3344b9ad928SBarry Smith Options Database Keys: 3354b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 3364b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 3374b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 3384b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 3394b9ad928SBarry Smith - -pc_sor_local_backward - Activates local backward version 3404b9ad928SBarry Smith 3414b9ad928SBarry Smith Notes: 3424b9ad928SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 3434b9ad928SBarry Smith which can be chosen with the option 3444b9ad928SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 3454b9ad928SBarry Smith 3464b9ad928SBarry Smith Level: intermediate 3474b9ad928SBarry Smith 3484b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric 3494b9ad928SBarry Smith 3504b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega() 3514b9ad928SBarry Smith @*/ 3527087cfbeSBarry Smith PetscErrorCode PCSORSetSymmetric(PC pc,MatSORType flag) 3534b9ad928SBarry Smith { 3544ac538c5SBarry Smith PetscErrorCode ierr; 3554b9ad928SBarry Smith 3564b9ad928SBarry Smith PetscFunctionBegin; 3570700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 358c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,flag,2); 3594ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSORSetSymmetric_C",(PC,MatSORType),(pc,flag));CHKERRQ(ierr); 3604b9ad928SBarry Smith PetscFunctionReturn(0); 3614b9ad928SBarry Smith } 3624b9ad928SBarry Smith 3634b9ad928SBarry Smith #undef __FUNCT__ 3644b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega" 3654b9ad928SBarry Smith /*@ 3664b9ad928SBarry Smith PCSORSetOmega - Sets the SOR relaxation coefficient, omega 3674b9ad928SBarry Smith (where omega = 1.0 by default). 3684b9ad928SBarry Smith 3693f9fe445SBarry Smith Logically Collective on PC 3704b9ad928SBarry Smith 3714b9ad928SBarry Smith Input Parameters: 3724b9ad928SBarry Smith + pc - the preconditioner context 3734b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2). 3744b9ad928SBarry Smith 3754b9ad928SBarry Smith Options Database Key: 3764b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 3774b9ad928SBarry Smith 3784b9ad928SBarry Smith Level: intermediate 3794b9ad928SBarry Smith 3804b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega 3814b9ad928SBarry Smith 3824b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega() 3834b9ad928SBarry Smith @*/ 3847087cfbeSBarry Smith PetscErrorCode PCSORSetOmega(PC pc,PetscReal omega) 3854b9ad928SBarry Smith { 3864ac538c5SBarry Smith PetscErrorCode ierr; 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith PetscFunctionBegin; 389c5eb9154SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 390c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc,omega,2); 3914ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSORSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr); 3924b9ad928SBarry Smith PetscFunctionReturn(0); 3934b9ad928SBarry Smith } 3944b9ad928SBarry Smith 3954b9ad928SBarry Smith #undef __FUNCT__ 3964b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations" 3974b9ad928SBarry Smith /*@ 3984b9ad928SBarry Smith PCSORSetIterations - Sets the number of inner iterations to 3994b9ad928SBarry Smith be used by the SOR preconditioner. The default is 1. 4004b9ad928SBarry Smith 4013f9fe445SBarry Smith Logically Collective on PC 4024b9ad928SBarry Smith 4034b9ad928SBarry Smith Input Parameters: 4044b9ad928SBarry Smith + pc - the preconditioner context 4054b9ad928SBarry Smith . lits - number of local iterations, smoothings over just variables on processor 4064b9ad928SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 4074b9ad928SBarry Smith 4084b9ad928SBarry Smith Options Database Key: 4094b9ad928SBarry Smith + -pc_sor_its <its> - Sets number of iterations 4104b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith Level: intermediate 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith Notes: When run on one processor the number of smoothings is lits*its 4154b9ad928SBarry Smith 4164b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations 4174b9ad928SBarry Smith 4184b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric() 4194b9ad928SBarry Smith @*/ 4207087cfbeSBarry Smith PetscErrorCode PCSORSetIterations(PC pc,PetscInt its,PetscInt lits) 4214b9ad928SBarry Smith { 4224ac538c5SBarry Smith PetscErrorCode ierr; 4234b9ad928SBarry Smith 4244b9ad928SBarry Smith PetscFunctionBegin; 4250700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 426c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,its,2); 4274ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCSORSetIterations_C",(PC,PetscInt,PetscInt),(pc,its,lits));CHKERRQ(ierr); 4284b9ad928SBarry Smith PetscFunctionReturn(0); 4294b9ad928SBarry Smith } 4304b9ad928SBarry Smith 4314b9ad928SBarry Smith /*MC 4324b9ad928SBarry Smith PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning 4334b9ad928SBarry Smith 4344b9ad928SBarry Smith Options Database Keys: 4354b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 4364b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 437a9510f2eSBarry Smith . -pc_sor_forward - Activates forward version 4384b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 439a9510f2eSBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version (default version) 4404b9ad928SBarry Smith . -pc_sor_local_backward - Activates local backward version 4414b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 442422a814eSBarry Smith . -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal 443a9510f2eSBarry Smith . -pc_sor_its <its> - Sets number of iterations (default 1) 444a9510f2eSBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations (default 1) 4454b9ad928SBarry Smith 4464b9ad928SBarry Smith Level: beginner 4474b9ad928SBarry Smith 4484b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel 4494b9ad928SBarry Smith 45037a17b4dSBarry Smith Notes: Only implemented for the AIJ and SeqBAIJ matrix formats. 4514b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 4524b9ad928SBarry Smith Jacobi with SOR on each block. 4534b9ad928SBarry Smith 454422a814eSBarry Smith For AIJ matrix if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that 455422a814eSBarry Smith zero will be used and hence the KSPSolve() will terminate with KSP_DIVERGED_NANORIF. If the option 456422a814eSBarry Smith KSPSetErrorIfNotConverged() or -ksp_error_if_not_converged the code will terminate as soon as it detects the 457422a814eSBarry Smith zero pivot. 458422a814eSBarry Smith 45937a17b4dSBarry Smith For SeqBAIJ matrices this implements point-block SOR, but the omega, its, lits options are not supported. 46037a17b4dSBarry Smith 461422a814eSBarry Smith For SeqBAIJ the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected 462422a814eSBarry Smith the computation is stopped with an error 463422a814eSBarry Smith 4644b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 4654b9ad928SBarry Smith PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT 4664b9ad928SBarry Smith M*/ 4674b9ad928SBarry Smith 4684b9ad928SBarry Smith #undef __FUNCT__ 4694b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR" 4708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc) 4714b9ad928SBarry Smith { 472dfbe8321SBarry Smith PetscErrorCode ierr; 4734b9ad928SBarry Smith PC_SOR *jac; 4744b9ad928SBarry Smith 4754b9ad928SBarry Smith PetscFunctionBegin; 476b00a9115SJed Brown ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 4774b9ad928SBarry Smith 4784b9ad928SBarry Smith pc->ops->apply = PCApply_SOR; 4794b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_SOR; 4804b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SOR; 4814b9ad928SBarry Smith pc->ops->setup = 0; 4824b9ad928SBarry Smith pc->ops->view = PCView_SOR; 4834b9ad928SBarry Smith pc->ops->destroy = PCDestroy_SOR; 4844b9ad928SBarry Smith pc->data = (void*)jac; 485d9bc8e36SBarry Smith jac->sym = SOR_LOCAL_SYMMETRIC_SWEEP; 4864b9ad928SBarry Smith jac->omega = 1.0; 48796fc60bcSBarry Smith jac->fshift = 0.0; 4884b9ad928SBarry Smith jac->its = 1; 4894b9ad928SBarry Smith jac->lits = 1; 4904b9ad928SBarry Smith 491bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetSymmetric_C",PCSORSetSymmetric_SOR);CHKERRQ(ierr); 492bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetOmega_C",PCSORSetOmega_SOR);CHKERRQ(ierr); 493bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORSetIterations_C",PCSORSetIterations_SOR);CHKERRQ(ierr); 494c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetSymmetric_C",PCSORGetSymmetric_SOR);CHKERRQ(ierr); 495c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetOmega_C",PCSORGetOmega_SOR);CHKERRQ(ierr); 496c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCSORGetIterations_C",PCSORGetIterations_SOR);CHKERRQ(ierr); 4974b9ad928SBarry Smith PetscFunctionReturn(0); 4984b9ad928SBarry Smith } 4994b9ad928SBarry Smith 5004b9ad928SBarry Smith 5014b9ad928SBarry Smith 5024b9ad928SBarry Smith 5034b9ad928SBarry Smith 504