14b9ad928SBarry Smith /* 24b9ad928SBarry Smith Defines a (S)SOR preconditioner for any Mat implementation 34b9ad928SBarry Smith */ 4af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 54b9ad928SBarry Smith 64b9ad928SBarry Smith typedef struct { 7c1ac3661SBarry Smith PetscInt its; /* inner iterations, number of sweeps */ 8c1ac3661SBarry Smith PetscInt lits; /* local inner iterations, number of sweeps applied by the local matrix mat->A */ 94b9ad928SBarry Smith MatSORType sym; /* forward, reverse, symmetric etc. */ 104b9ad928SBarry Smith PetscReal omega; 1129c1d7e0SHong Zhang PetscReal fshift; 124b9ad928SBarry Smith } PC_SOR; 134b9ad928SBarry Smith 14d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_SOR(PC pc) 15d71ae5a4SJacob Faibussowitsch { 164b9ad928SBarry Smith PetscFunctionBegin; 172e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetSymmetric_C", NULL)); 182e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetOmega_C", NULL)); 192e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetIterations_C", NULL)); 202e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetSymmetric_C", NULL)); 212e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetOmega_C", NULL)); 222e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetIterations_C", NULL)); 239566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 244b9ad928SBarry Smith PetscFunctionReturn(0); 254b9ad928SBarry Smith } 264b9ad928SBarry Smith 27d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_SOR(PC pc, Vec x, Vec y) 28d71ae5a4SJacob Faibussowitsch { 294b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 30c1ac3661SBarry Smith PetscInt flag = jac->sym | SOR_ZERO_INITIAL_GUESS; 314b9ad928SBarry Smith 324b9ad928SBarry Smith PetscFunctionBegin; 339566063dSJacob Faibussowitsch PetscCall(MatSOR(pc->pmat, x, jac->omega, (MatSORType)flag, jac->fshift, jac->its, jac->lits, y)); 349566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, (MatFactorError *)&pc->failedreason)); 354b9ad928SBarry Smith PetscFunctionReturn(0); 364b9ad928SBarry Smith } 374b9ad928SBarry Smith 38d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_SOR(PC pc, Vec x, Vec y) 39d71ae5a4SJacob Faibussowitsch { 409d2471e0SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 419d2471e0SBarry Smith PetscInt flag = jac->sym | SOR_ZERO_INITIAL_GUESS; 429d2471e0SBarry Smith PetscBool set, sym; 439d2471e0SBarry Smith 449d2471e0SBarry Smith PetscFunctionBegin; 459566063dSJacob Faibussowitsch PetscCall(MatIsSymmetricKnown(pc->pmat, &set, &sym)); 467827d75bSBarry Smith PetscCheck(set && sym && (jac->sym == SOR_SYMMETRIC_SWEEP || jac->sym == SOR_LOCAL_SYMMETRIC_SWEEP), PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Can only apply transpose of SOR if matrix is symmetric and sweep is symmetric"); 479566063dSJacob Faibussowitsch PetscCall(MatSOR(pc->pmat, x, jac->omega, (MatSORType)flag, jac->fshift, jac->its, jac->lits, y)); 489566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, (MatFactorError *)&pc->failedreason)); 499d2471e0SBarry Smith PetscFunctionReturn(0); 509d2471e0SBarry Smith } 519d2471e0SBarry Smith 52d71ae5a4SJacob Faibussowitsch 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) 53d71ae5a4SJacob Faibussowitsch { 544b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 557319c654SBarry Smith MatSORType stype = jac->sym; 564b9ad928SBarry Smith 574b9ad928SBarry Smith PetscFunctionBegin; 5835cb6cd3SPierre Jolivet PetscCall(PetscInfo(pc, "Warning, convergence criteria ignored, using %" PetscInt_FMT " iterations\n", its)); 592fa5cd67SKarl Rupp if (guesszero) stype = (MatSORType)(stype | SOR_ZERO_INITIAL_GUESS); 609566063dSJacob Faibussowitsch PetscCall(MatSOR(pc->pmat, b, jac->omega, stype, jac->fshift, its * jac->its, jac->lits, y)); 619566063dSJacob Faibussowitsch PetscCall(MatFactorGetError(pc->pmat, (MatFactorError *)&pc->failedreason)); 624d0a8057SBarry Smith *outits = its; 634d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ITS; 644b9ad928SBarry Smith PetscFunctionReturn(0); 654b9ad928SBarry Smith } 664b9ad928SBarry Smith 67d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSetFromOptions_SOR(PC pc, PetscOptionItems *PetscOptionsObject) 68d71ae5a4SJacob Faibussowitsch { 694b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 70ace3abfcSBarry Smith PetscBool flg; 714b9ad928SBarry Smith 724b9ad928SBarry Smith PetscFunctionBegin; 73d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "(S)SOR options"); 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_sor_omega", "relaxation factor (0 < omega < 2)", "PCSORSetOmega", jac->omega, &jac->omega, NULL)); 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_sor_diagonal_shift", "Add to the diagonal entries", "", jac->fshift, &jac->fshift, NULL)); 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_sor_its", "number of inner SOR iterations", "PCSORSetIterations", jac->its, &jac->its, NULL)); 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_sor_lits", "number of local inner SOR iterations", "PCSORSetIterations", jac->lits, &jac->lits, NULL)); 789566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroupBegin("-pc_sor_symmetric", "SSOR, not SOR", "PCSORSetSymmetric", &flg)); 799566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_SYMMETRIC_SWEEP)); 809566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_backward", "use backward sweep instead of forward", "PCSORSetSymmetric", &flg)); 819566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_BACKWARD_SWEEP)); 829566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_forward", "use forward sweep", "PCSORSetSymmetric", &flg)); 839566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_FORWARD_SWEEP)); 849566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_local_symmetric", "use SSOR separately on each processor", "PCSORSetSymmetric", &flg)); 859566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_LOCAL_SYMMETRIC_SWEEP)); 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_local_backward", "use backward sweep locally", "PCSORSetSymmetric", &flg)); 879566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_LOCAL_BACKWARD_SWEEP)); 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroupEnd("-pc_sor_local_forward", "use forward sweep locally", "PCSORSetSymmetric", &flg)); 899566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_LOCAL_FORWARD_SWEEP)); 90d0609cedSBarry Smith PetscOptionsHeadEnd(); 914b9ad928SBarry Smith PetscFunctionReturn(0); 924b9ad928SBarry Smith } 934b9ad928SBarry Smith 94d71ae5a4SJacob Faibussowitsch PetscErrorCode PCView_SOR(PC pc, PetscViewer viewer) 95d71ae5a4SJacob Faibussowitsch { 964b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 974b9ad928SBarry Smith MatSORType sym = jac->sym; 982fc52814SBarry Smith const char *sortype; 99ace3abfcSBarry Smith PetscBool iascii; 1004b9ad928SBarry Smith 1014b9ad928SBarry Smith PetscFunctionBegin; 1029566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 10332077d6dSBarry Smith if (iascii) { 1049566063dSJacob Faibussowitsch if (sym & SOR_ZERO_INITIAL_GUESS) PetscCall(PetscViewerASCIIPrintf(viewer, " zero initial guess\n")); 1054b9ad928SBarry Smith if (sym == SOR_APPLY_UPPER) sortype = "apply_upper"; 1064b9ad928SBarry Smith else if (sym == SOR_APPLY_LOWER) sortype = "apply_lower"; 1074b9ad928SBarry Smith else if (sym & SOR_EISENSTAT) sortype = "Eisenstat"; 108db4deed7SKarl Rupp else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP) sortype = "symmetric"; 1094b9ad928SBarry Smith else if (sym & SOR_BACKWARD_SWEEP) sortype = "backward"; 1104b9ad928SBarry Smith else if (sym & SOR_FORWARD_SWEEP) sortype = "forward"; 111db4deed7SKarl Rupp else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) 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"; 11563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type = %s, iterations = %" PetscInt_FMT ", local iterations = %" PetscInt_FMT ", omega = %g\n", sortype, jac->its, jac->lits, (double)jac->omega)); 1164b9ad928SBarry Smith } 1174b9ad928SBarry Smith PetscFunctionReturn(0); 1184b9ad928SBarry Smith } 1194b9ad928SBarry Smith 120d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORSetSymmetric_SOR(PC pc, MatSORType flag) 121d71ae5a4SJacob Faibussowitsch { 122c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 1234b9ad928SBarry Smith 1244b9ad928SBarry Smith PetscFunctionBegin; 1254b9ad928SBarry Smith jac->sym = flag; 1264b9ad928SBarry Smith PetscFunctionReturn(0); 1274b9ad928SBarry Smith } 1284b9ad928SBarry Smith 129d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORSetOmega_SOR(PC pc, PetscReal omega) 130d71ae5a4SJacob Faibussowitsch { 131c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 1324b9ad928SBarry Smith 1334b9ad928SBarry Smith PetscFunctionBegin; 1342472a847SBarry Smith PetscCheck(omega > 0.0 && omega < 2.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Relaxation out of range"); 1354b9ad928SBarry Smith jac->omega = omega; 1364b9ad928SBarry Smith PetscFunctionReturn(0); 1374b9ad928SBarry Smith } 1384b9ad928SBarry Smith 139d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORSetIterations_SOR(PC pc, PetscInt its, PetscInt lits) 140d71ae5a4SJacob Faibussowitsch { 141c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 1424b9ad928SBarry Smith 1434b9ad928SBarry Smith PetscFunctionBegin; 1444b9ad928SBarry Smith jac->its = its; 1454b9ad928SBarry Smith jac->lits = lits; 1464b9ad928SBarry Smith PetscFunctionReturn(0); 1474b9ad928SBarry Smith } 1484b9ad928SBarry Smith 149d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORGetSymmetric_SOR(PC pc, MatSORType *flag) 150d71ae5a4SJacob Faibussowitsch { 151c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 152c60c7ad4SBarry Smith 153c60c7ad4SBarry Smith PetscFunctionBegin; 154c60c7ad4SBarry Smith *flag = jac->sym; 155c60c7ad4SBarry Smith PetscFunctionReturn(0); 156c60c7ad4SBarry Smith } 157c60c7ad4SBarry Smith 158d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORGetOmega_SOR(PC pc, PetscReal *omega) 159d71ae5a4SJacob Faibussowitsch { 160c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 161c60c7ad4SBarry Smith 162c60c7ad4SBarry Smith PetscFunctionBegin; 163c60c7ad4SBarry Smith *omega = jac->omega; 164c60c7ad4SBarry Smith PetscFunctionReturn(0); 165c60c7ad4SBarry Smith } 166c60c7ad4SBarry Smith 167d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORGetIterations_SOR(PC pc, PetscInt *its, PetscInt *lits) 168d71ae5a4SJacob Faibussowitsch { 169c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 170c60c7ad4SBarry Smith 171c60c7ad4SBarry Smith PetscFunctionBegin; 172c60c7ad4SBarry Smith if (its) *its = jac->its; 173c60c7ad4SBarry Smith if (lits) *lits = jac->lits; 174c60c7ad4SBarry Smith PetscFunctionReturn(0); 175c60c7ad4SBarry Smith } 176c60c7ad4SBarry Smith 177c60c7ad4SBarry Smith /*@ 1781ff2113eSBarry Smith PCSORGetSymmetric - Gets the form the SOR preconditioner is using; backward, or forward relaxation. The local variants perform SOR on 179c60c7ad4SBarry Smith each processor. By default forward relaxation is used. 180c60c7ad4SBarry Smith 181c3339decSBarry Smith Logically Collective 182c60c7ad4SBarry Smith 183c60c7ad4SBarry Smith Input Parameter: 184c60c7ad4SBarry Smith . pc - the preconditioner context 185c60c7ad4SBarry Smith 186c60c7ad4SBarry Smith Output Parameter: 187c60c7ad4SBarry Smith . flag - one of the following 188c60c7ad4SBarry Smith .vb 189c60c7ad4SBarry Smith SOR_FORWARD_SWEEP 190c60c7ad4SBarry Smith SOR_BACKWARD_SWEEP 191c60c7ad4SBarry Smith SOR_SYMMETRIC_SWEEP 192c60c7ad4SBarry Smith SOR_LOCAL_FORWARD_SWEEP 193c60c7ad4SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 194c60c7ad4SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 195c60c7ad4SBarry Smith .ve 196c60c7ad4SBarry Smith 197c60c7ad4SBarry Smith Options Database Keys: 198c60c7ad4SBarry Smith + -pc_sor_symmetric - Activates symmetric version 199c60c7ad4SBarry Smith . -pc_sor_backward - Activates backward version 200c60c7ad4SBarry Smith . -pc_sor_local_forward - Activates local forward version 201c60c7ad4SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 202c60c7ad4SBarry Smith - -pc_sor_local_backward - Activates local backward version 203c60c7ad4SBarry Smith 204f1580f4eSBarry Smith Note: 205f1580f4eSBarry Smith To use the Eisenstat trick with SSOR, employ the `PCEISENSTAT` preconditioner, 206c60c7ad4SBarry Smith which can be chosen with the option 207c60c7ad4SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 208c60c7ad4SBarry Smith 209c60c7ad4SBarry Smith Level: intermediate 210c60c7ad4SBarry Smith 211f1580f4eSBarry Smith .seealso: `PCSOR`, `PCEisenstatSetOmega()`, `PCSORSetIterations()`, `PCSORSetOmega()`, `PCSORSetSymmetric()` 212c60c7ad4SBarry Smith @*/ 213d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORGetSymmetric(PC pc, MatSORType *flag) 214d71ae5a4SJacob Faibussowitsch { 215c60c7ad4SBarry Smith PetscFunctionBegin; 216c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 217cac4c232SBarry Smith PetscUseMethod(pc, "PCSORGetSymmetric_C", (PC, MatSORType *), (pc, flag)); 218c60c7ad4SBarry Smith PetscFunctionReturn(0); 219c60c7ad4SBarry Smith } 220c60c7ad4SBarry Smith 221c60c7ad4SBarry Smith /*@ 222c60c7ad4SBarry Smith PCSORGetOmega - Gets the SOR relaxation coefficient, omega 223c60c7ad4SBarry Smith (where omega = 1.0 by default). 224c60c7ad4SBarry Smith 225c3339decSBarry Smith Logically Collective 226c60c7ad4SBarry Smith 227c60c7ad4SBarry Smith Input Parameter: 228c60c7ad4SBarry Smith . pc - the preconditioner context 229c60c7ad4SBarry Smith 230c60c7ad4SBarry Smith Output Parameter: 231c60c7ad4SBarry Smith . omega - relaxation coefficient (0 < omega < 2). 232c60c7ad4SBarry Smith 233c60c7ad4SBarry Smith Options Database Key: 234c60c7ad4SBarry Smith . -pc_sor_omega <omega> - Sets omega 235c60c7ad4SBarry Smith 236c60c7ad4SBarry Smith Level: intermediate 237c60c7ad4SBarry Smith 238f1580f4eSBarry Smith .seealso: `PCSOR`, `PCSORSetSymmetric()`, `PCSORSetIterations()`, `PCEisenstatSetOmega()`, `PCSORSetOmega()` 239c60c7ad4SBarry Smith @*/ 240d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORGetOmega(PC pc, PetscReal *omega) 241d71ae5a4SJacob Faibussowitsch { 242c60c7ad4SBarry Smith PetscFunctionBegin; 243c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 244cac4c232SBarry Smith PetscUseMethod(pc, "PCSORGetOmega_C", (PC, PetscReal *), (pc, omega)); 245c60c7ad4SBarry Smith PetscFunctionReturn(0); 246c60c7ad4SBarry Smith } 247c60c7ad4SBarry Smith 248c60c7ad4SBarry Smith /*@ 249c60c7ad4SBarry Smith PCSORGetIterations - Gets the number of inner iterations to 250c60c7ad4SBarry Smith be used by the SOR preconditioner. The default is 1. 251c60c7ad4SBarry Smith 252c3339decSBarry Smith Logically Collective 253c60c7ad4SBarry Smith 254c60c7ad4SBarry Smith Input Parameter: 255c60c7ad4SBarry Smith . pc - the preconditioner context 256c60c7ad4SBarry Smith 257d8d19677SJose E. Roman Output Parameters: 258c60c7ad4SBarry Smith + lits - number of local iterations, smoothings over just variables on processor 259c60c7ad4SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 260c60c7ad4SBarry Smith 261f1580f4eSBarry Smith Options Database Keys: 262c60c7ad4SBarry Smith + -pc_sor_its <its> - Sets number of iterations 263c60c7ad4SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 264c60c7ad4SBarry Smith 265c60c7ad4SBarry Smith Level: intermediate 266c60c7ad4SBarry Smith 267f1580f4eSBarry Smith Note: 26895452b02SPatrick Sanan When run on one processor the number of smoothings is lits*its 269c60c7ad4SBarry Smith 270f1580f4eSBarry Smith .seealso: `PCSOR`, `PCSORSetOmega()`, `PCSORSetSymmetric()`, `PCSORSetIterations()` 271c60c7ad4SBarry Smith @*/ 272d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORGetIterations(PC pc, PetscInt *its, PetscInt *lits) 273d71ae5a4SJacob Faibussowitsch { 274c60c7ad4SBarry Smith PetscFunctionBegin; 275c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 276cac4c232SBarry Smith PetscUseMethod(pc, "PCSORGetIterations_C", (PC, PetscInt *, PetscInt *), (pc, its, lits)); 277c60c7ad4SBarry Smith PetscFunctionReturn(0); 278c60c7ad4SBarry Smith } 279c60c7ad4SBarry Smith 2804b9ad928SBarry Smith /*@ 2814b9ad928SBarry Smith PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR), 2824b9ad928SBarry Smith backward, or forward relaxation. The local variants perform SOR on 2834b9ad928SBarry Smith each processor. By default forward relaxation is used. 2844b9ad928SBarry Smith 285c3339decSBarry Smith Logically Collective 2864b9ad928SBarry Smith 2874b9ad928SBarry Smith Input Parameters: 2884b9ad928SBarry Smith + pc - the preconditioner context 2894b9ad928SBarry Smith - flag - one of the following 2904b9ad928SBarry Smith .vb 2914b9ad928SBarry Smith SOR_FORWARD_SWEEP 2924b9ad928SBarry Smith SOR_BACKWARD_SWEEP 2934b9ad928SBarry Smith SOR_SYMMETRIC_SWEEP 2944b9ad928SBarry Smith SOR_LOCAL_FORWARD_SWEEP 2954b9ad928SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 2964b9ad928SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 2974b9ad928SBarry Smith .ve 2984b9ad928SBarry Smith 2994b9ad928SBarry Smith Options Database Keys: 3004b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 3014b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 3024b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 3034b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 3044b9ad928SBarry Smith - -pc_sor_local_backward - Activates local backward version 3054b9ad928SBarry Smith 306f1580f4eSBarry Smith Note: 3074b9ad928SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 3084b9ad928SBarry Smith which can be chosen with the option 3094b9ad928SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 3104b9ad928SBarry Smith 3114b9ad928SBarry Smith Level: intermediate 3124b9ad928SBarry Smith 313f1580f4eSBarry Smith .seealso: `PCSOR`, `PCEisenstatSetOmega()`, `PCSORSetIterations()`, `PCSORSetOmega()` 3144b9ad928SBarry Smith @*/ 315d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORSetSymmetric(PC pc, MatSORType flag) 316d71ae5a4SJacob Faibussowitsch { 3174b9ad928SBarry Smith PetscFunctionBegin; 3180700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 319c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, flag, 2); 320cac4c232SBarry Smith PetscTryMethod(pc, "PCSORSetSymmetric_C", (PC, MatSORType), (pc, flag)); 3214b9ad928SBarry Smith PetscFunctionReturn(0); 3224b9ad928SBarry Smith } 3234b9ad928SBarry Smith 3244b9ad928SBarry Smith /*@ 3254b9ad928SBarry Smith PCSORSetOmega - Sets the SOR relaxation coefficient, omega 3264b9ad928SBarry Smith (where omega = 1.0 by default). 3274b9ad928SBarry Smith 328c3339decSBarry Smith Logically Collective 3294b9ad928SBarry Smith 3304b9ad928SBarry Smith Input Parameters: 3314b9ad928SBarry Smith + pc - the preconditioner context 3324b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2). 3334b9ad928SBarry Smith 3344b9ad928SBarry Smith Options Database Key: 3354b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 3364b9ad928SBarry Smith 3374b9ad928SBarry Smith Level: intermediate 3384b9ad928SBarry Smith 339485168efSMatthew Knepley Note: 340f1580f4eSBarry Smith If omega != 1, you will need to set the `MAT_USE_INODE`S option to `PETSC_FALSE` on the matrix. 341485168efSMatthew Knepley 342f1580f4eSBarry Smith .seealso: `PCSOR`, `PCSORSetSymmetric()`, `PCSORSetIterations()`, `PCEisenstatSetOmega()`, `MatSetOption()` 3434b9ad928SBarry Smith @*/ 344d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORSetOmega(PC pc, PetscReal omega) 345d71ae5a4SJacob Faibussowitsch { 3464b9ad928SBarry Smith PetscFunctionBegin; 347c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 348c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, omega, 2); 349cac4c232SBarry Smith PetscTryMethod(pc, "PCSORSetOmega_C", (PC, PetscReal), (pc, omega)); 3504b9ad928SBarry Smith PetscFunctionReturn(0); 3514b9ad928SBarry Smith } 3524b9ad928SBarry Smith 3534b9ad928SBarry Smith /*@ 3544b9ad928SBarry Smith PCSORSetIterations - Sets the number of inner iterations to 3554b9ad928SBarry Smith be used by the SOR preconditioner. The default is 1. 3564b9ad928SBarry Smith 357c3339decSBarry Smith Logically Collective 3584b9ad928SBarry Smith 3594b9ad928SBarry Smith Input Parameters: 3604b9ad928SBarry Smith + pc - the preconditioner context 3614b9ad928SBarry Smith . lits - number of local iterations, smoothings over just variables on processor 3624b9ad928SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 3634b9ad928SBarry Smith 364f1580f4eSBarry Smith Options Database Keys: 3654b9ad928SBarry Smith + -pc_sor_its <its> - Sets number of iterations 3664b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 3674b9ad928SBarry Smith 3684b9ad928SBarry Smith Level: intermediate 3694b9ad928SBarry Smith 370f1580f4eSBarry Smith Note: 37195452b02SPatrick Sanan When run on one processor the number of smoothings is lits*its 3724b9ad928SBarry Smith 373f1580f4eSBarry Smith .seealso: `PCSOR`, `PCSORSetOmega()`, `PCSORSetSymmetric()` 3744b9ad928SBarry Smith @*/ 375d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORSetIterations(PC pc, PetscInt its, PetscInt lits) 376d71ae5a4SJacob Faibussowitsch { 3774b9ad928SBarry Smith PetscFunctionBegin; 3780700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 379c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, its, 2); 380cac4c232SBarry Smith PetscTryMethod(pc, "PCSORSetIterations_C", (PC, PetscInt, PetscInt), (pc, its, lits)); 3814b9ad928SBarry Smith PetscFunctionReturn(0); 3824b9ad928SBarry Smith } 3834b9ad928SBarry Smith 3844b9ad928SBarry Smith /*MC 3854b9ad928SBarry Smith PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning 3864b9ad928SBarry Smith 3874b9ad928SBarry Smith Options Database Keys: 3884b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 3894b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 390a9510f2eSBarry Smith . -pc_sor_forward - Activates forward version 3914b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 392a9510f2eSBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version (default version) 3934b9ad928SBarry Smith . -pc_sor_local_backward - Activates local backward version 3944b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 395422a814eSBarry Smith . -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal 396a9510f2eSBarry Smith . -pc_sor_its <its> - Sets number of iterations (default 1) 397a9510f2eSBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations (default 1) 3984b9ad928SBarry Smith 3994b9ad928SBarry Smith Level: beginner 4004b9ad928SBarry Smith 40195452b02SPatrick Sanan Notes: 402f1580f4eSBarry Smith Only implemented for the `MATAIJ` and `MATSEQBAIJ` matrix formats. 403f1580f4eSBarry Smith 4044b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 4054b9ad928SBarry Smith Jacobi with SOR on each block. 4064b9ad928SBarry Smith 407f1580f4eSBarry Smith For `MATAIJ` matrices if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that 408*48773899SPierre Jolivet zero will be used and hence the `KSPSolve()` will terminate with `KSP_DIVERGED_NANORINF`. If the option 409f1580f4eSBarry Smith `KSPSetErrorIfNotConverged()` or -ksp_error_if_not_converged the code will terminate as soon as it detects the 410422a814eSBarry Smith zero pivot. 411422a814eSBarry Smith 412f1580f4eSBarry Smith For `MATSEQBAIJ` matrices this implements point-block SOR, but the omega, its, lits options are not supported. 41337a17b4dSBarry Smith 414f1580f4eSBarry Smith For `MATSEQBAIJ` the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected 415422a814eSBarry Smith the computation is stopped with an error 416422a814eSBarry Smith 417f1580f4eSBarry Smith If used with `KSPRICHARDSON` and no monitors the convergence test is skipped to improve speed, thus it always iterates 418f1580f4eSBarry Smith the maximum number of iterations you've selected for `KSP`. It is usually used in this mode as a smoother for multigrid. 41967991ba4SBarry Smith 420f1580f4eSBarry Smith If omega != 1, you will need to set the `MAT_USE_INODES` option to `PETSC_FALSE` on the matrix. 421485168efSMatthew Knepley 422f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCJACOBI`, 423db781477SPatrick Sanan `PCSORSetIterations()`, `PCSORSetSymmetric()`, `PCSORSetOmega()`, `PCEISENSTAT`, `MatSetOption()` 4244b9ad928SBarry Smith M*/ 4254b9ad928SBarry Smith 426d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc) 427d71ae5a4SJacob Faibussowitsch { 4284b9ad928SBarry Smith PC_SOR *jac; 4294b9ad928SBarry Smith 4304b9ad928SBarry Smith PetscFunctionBegin; 4314dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 4324b9ad928SBarry Smith 4334b9ad928SBarry Smith pc->ops->apply = PCApply_SOR; 4349d2471e0SBarry Smith pc->ops->applytranspose = PCApplyTranspose_SOR; 4354b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_SOR; 4364b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SOR; 4370a545947SLisandro Dalcin pc->ops->setup = NULL; 4384b9ad928SBarry Smith pc->ops->view = PCView_SOR; 4394b9ad928SBarry Smith pc->ops->destroy = PCDestroy_SOR; 4404b9ad928SBarry Smith pc->data = (void *)jac; 441d9bc8e36SBarry Smith jac->sym = SOR_LOCAL_SYMMETRIC_SWEEP; 4424b9ad928SBarry Smith jac->omega = 1.0; 44396fc60bcSBarry Smith jac->fshift = 0.0; 4444b9ad928SBarry Smith jac->its = 1; 4454b9ad928SBarry Smith jac->lits = 1; 4464b9ad928SBarry Smith 4479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetSymmetric_C", PCSORSetSymmetric_SOR)); 4489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetOmega_C", PCSORSetOmega_SOR)); 4499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetIterations_C", PCSORSetIterations_SOR)); 4509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetSymmetric_C", PCSORGetSymmetric_SOR)); 4519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetOmega_C", PCSORGetOmega_SOR)); 4529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetIterations_C", PCSORGetIterations_SOR)); 4534b9ad928SBarry Smith PetscFunctionReturn(0); 4544b9ad928SBarry Smith } 455