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)); 243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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)); 353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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)); 493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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; 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 654b9ad928SBarry Smith } 664b9ad928SBarry Smith 6766976f2fSJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_SOR(PC pc, PetscOptionItems *PetscOptionsObject) 68d71ae5a4SJacob Faibussowitsch { 694b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 70ace3abfcSBarry Smith PetscBool flg; 71d8e4f26eSJose E. Roman PetscReal omega; 724b9ad928SBarry Smith 734b9ad928SBarry Smith PetscFunctionBegin; 74d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "(S)SOR options"); 75d8e4f26eSJose E. Roman PetscCall(PetscOptionsReal("-pc_sor_omega", "relaxation factor (0 < omega < 2)", "PCSORSetOmega", jac->omega, &omega, &flg)); 76d8e4f26eSJose E. Roman if (flg) PetscCall(PCSORSetOmega(pc, omega)); 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_sor_diagonal_shift", "Add to the diagonal entries", "", jac->fshift, &jac->fshift, NULL)); 789566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_sor_its", "number of inner SOR iterations", "PCSORSetIterations", jac->its, &jac->its, NULL)); 799566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_sor_lits", "number of local inner SOR iterations", "PCSORSetIterations", jac->lits, &jac->lits, NULL)); 809566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroupBegin("-pc_sor_symmetric", "SSOR, not SOR", "PCSORSetSymmetric", &flg)); 819566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_SYMMETRIC_SWEEP)); 829566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_backward", "use backward sweep instead of forward", "PCSORSetSymmetric", &flg)); 839566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_BACKWARD_SWEEP)); 849566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_forward", "use forward sweep", "PCSORSetSymmetric", &flg)); 859566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_FORWARD_SWEEP)); 869566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_local_symmetric", "use SSOR separately on each processor", "PCSORSetSymmetric", &flg)); 879566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_LOCAL_SYMMETRIC_SWEEP)); 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroup("-pc_sor_local_backward", "use backward sweep locally", "PCSORSetSymmetric", &flg)); 899566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_LOCAL_BACKWARD_SWEEP)); 909566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoolGroupEnd("-pc_sor_local_forward", "use forward sweep locally", "PCSORSetSymmetric", &flg)); 919566063dSJacob Faibussowitsch if (flg) PetscCall(PCSORSetSymmetric(pc, SOR_LOCAL_FORWARD_SWEEP)); 92d0609cedSBarry Smith PetscOptionsHeadEnd(); 933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 944b9ad928SBarry Smith } 954b9ad928SBarry Smith 9666976f2fSJacob Faibussowitsch static PetscErrorCode PCView_SOR(PC pc, PetscViewer viewer) 97d71ae5a4SJacob Faibussowitsch { 984b9ad928SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 994b9ad928SBarry Smith MatSORType sym = jac->sym; 1002fc52814SBarry Smith const char *sortype; 101ace3abfcSBarry Smith PetscBool iascii; 1024b9ad928SBarry Smith 1034b9ad928SBarry Smith PetscFunctionBegin; 1049566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 10532077d6dSBarry Smith if (iascii) { 1069566063dSJacob Faibussowitsch if (sym & SOR_ZERO_INITIAL_GUESS) PetscCall(PetscViewerASCIIPrintf(viewer, " zero initial guess\n")); 1074b9ad928SBarry Smith if (sym == SOR_APPLY_UPPER) sortype = "apply_upper"; 1084b9ad928SBarry Smith else if (sym == SOR_APPLY_LOWER) sortype = "apply_lower"; 1094b9ad928SBarry Smith else if (sym & SOR_EISENSTAT) sortype = "Eisenstat"; 110db4deed7SKarl Rupp else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP) sortype = "symmetric"; 1114b9ad928SBarry Smith else if (sym & SOR_BACKWARD_SWEEP) sortype = "backward"; 1124b9ad928SBarry Smith else if (sym & SOR_FORWARD_SWEEP) sortype = "forward"; 113db4deed7SKarl Rupp else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) sortype = "local_symmetric"; 1144b9ad928SBarry Smith else if (sym & SOR_LOCAL_FORWARD_SWEEP) sortype = "local_forward"; 1154b9ad928SBarry Smith else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward"; 1164b9ad928SBarry Smith else sortype = "unknown"; 11763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type = %s, iterations = %" PetscInt_FMT ", local iterations = %" PetscInt_FMT ", omega = %g\n", sortype, jac->its, jac->lits, (double)jac->omega)); 1184b9ad928SBarry Smith } 1193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1204b9ad928SBarry Smith } 1214b9ad928SBarry Smith 122d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORSetSymmetric_SOR(PC pc, MatSORType flag) 123d71ae5a4SJacob Faibussowitsch { 124c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 1254b9ad928SBarry Smith 1264b9ad928SBarry Smith PetscFunctionBegin; 1274b9ad928SBarry Smith jac->sym = flag; 1283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1294b9ad928SBarry Smith } 1304b9ad928SBarry Smith 131d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORSetOmega_SOR(PC pc, PetscReal omega) 132d71ae5a4SJacob Faibussowitsch { 133c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 1344b9ad928SBarry Smith 1354b9ad928SBarry Smith PetscFunctionBegin; 1362472a847SBarry Smith PetscCheck(omega > 0.0 && omega < 2.0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Relaxation out of range"); 1374b9ad928SBarry Smith jac->omega = omega; 1383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1394b9ad928SBarry Smith } 1404b9ad928SBarry Smith 141d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORSetIterations_SOR(PC pc, PetscInt its, PetscInt lits) 142d71ae5a4SJacob Faibussowitsch { 143c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 1444b9ad928SBarry Smith 1454b9ad928SBarry Smith PetscFunctionBegin; 1464b9ad928SBarry Smith jac->its = its; 1474b9ad928SBarry Smith jac->lits = lits; 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1494b9ad928SBarry Smith } 1504b9ad928SBarry Smith 151d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORGetSymmetric_SOR(PC pc, MatSORType *flag) 152d71ae5a4SJacob Faibussowitsch { 153c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 154c60c7ad4SBarry Smith 155c60c7ad4SBarry Smith PetscFunctionBegin; 156c60c7ad4SBarry Smith *flag = jac->sym; 1573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 158c60c7ad4SBarry Smith } 159c60c7ad4SBarry Smith 160d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORGetOmega_SOR(PC pc, PetscReal *omega) 161d71ae5a4SJacob Faibussowitsch { 162c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 163c60c7ad4SBarry Smith 164c60c7ad4SBarry Smith PetscFunctionBegin; 165c60c7ad4SBarry Smith *omega = jac->omega; 1663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 167c60c7ad4SBarry Smith } 168c60c7ad4SBarry Smith 169d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSORGetIterations_SOR(PC pc, PetscInt *its, PetscInt *lits) 170d71ae5a4SJacob Faibussowitsch { 171c60c7ad4SBarry Smith PC_SOR *jac = (PC_SOR *)pc->data; 172c60c7ad4SBarry Smith 173c60c7ad4SBarry Smith PetscFunctionBegin; 174c60c7ad4SBarry Smith if (its) *its = jac->its; 175c60c7ad4SBarry Smith if (lits) *lits = jac->lits; 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 177c60c7ad4SBarry Smith } 178c60c7ad4SBarry Smith 179c60c7ad4SBarry Smith /*@ 1801ff2113eSBarry Smith PCSORGetSymmetric - Gets the form the SOR preconditioner is using; backward, or forward relaxation. The local variants perform SOR on 181c60c7ad4SBarry Smith each processor. By default forward relaxation is used. 182c60c7ad4SBarry Smith 183c3339decSBarry Smith Logically Collective 184c60c7ad4SBarry Smith 185c60c7ad4SBarry Smith Input Parameter: 186c60c7ad4SBarry Smith . pc - the preconditioner context 187c60c7ad4SBarry Smith 188c60c7ad4SBarry Smith Output Parameter: 189c60c7ad4SBarry Smith . flag - one of the following 190c60c7ad4SBarry Smith .vb 191c60c7ad4SBarry Smith SOR_FORWARD_SWEEP 192c60c7ad4SBarry Smith SOR_BACKWARD_SWEEP 193c60c7ad4SBarry Smith SOR_SYMMETRIC_SWEEP 194c60c7ad4SBarry Smith SOR_LOCAL_FORWARD_SWEEP 195c60c7ad4SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 196c60c7ad4SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 197c60c7ad4SBarry Smith .ve 198c60c7ad4SBarry Smith 199c60c7ad4SBarry Smith Options Database Keys: 200c60c7ad4SBarry Smith + -pc_sor_symmetric - Activates symmetric version 201c60c7ad4SBarry Smith . -pc_sor_backward - Activates backward version 202c60c7ad4SBarry Smith . -pc_sor_local_forward - Activates local forward version 203c60c7ad4SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 204c60c7ad4SBarry Smith - -pc_sor_local_backward - Activates local backward version 205c60c7ad4SBarry Smith 206f1580f4eSBarry Smith Note: 207f1580f4eSBarry Smith To use the Eisenstat trick with SSOR, employ the `PCEISENSTAT` preconditioner, 208c60c7ad4SBarry Smith which can be chosen with the option 209c60c7ad4SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 210c60c7ad4SBarry Smith 211c60c7ad4SBarry Smith Level: intermediate 212c60c7ad4SBarry Smith 213*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCEisenstatSetOmega()`, `PCSORSetIterations()`, `PCSORSetOmega()`, `PCSORSetSymmetric()` 214c60c7ad4SBarry Smith @*/ 215d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORGetSymmetric(PC pc, MatSORType *flag) 216d71ae5a4SJacob Faibussowitsch { 217c60c7ad4SBarry Smith PetscFunctionBegin; 218c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 219cac4c232SBarry Smith PetscUseMethod(pc, "PCSORGetSymmetric_C", (PC, MatSORType *), (pc, flag)); 2203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 221c60c7ad4SBarry Smith } 222c60c7ad4SBarry Smith 223c60c7ad4SBarry Smith /*@ 224c60c7ad4SBarry Smith PCSORGetOmega - Gets the SOR relaxation coefficient, omega 225c60c7ad4SBarry Smith (where omega = 1.0 by default). 226c60c7ad4SBarry Smith 227c3339decSBarry Smith Logically Collective 228c60c7ad4SBarry Smith 229c60c7ad4SBarry Smith Input Parameter: 230c60c7ad4SBarry Smith . pc - the preconditioner context 231c60c7ad4SBarry Smith 232c60c7ad4SBarry Smith Output Parameter: 233c60c7ad4SBarry Smith . omega - relaxation coefficient (0 < omega < 2). 234c60c7ad4SBarry Smith 235c60c7ad4SBarry Smith Options Database Key: 236c60c7ad4SBarry Smith . -pc_sor_omega <omega> - Sets omega 237c60c7ad4SBarry Smith 238c60c7ad4SBarry Smith Level: intermediate 239c60c7ad4SBarry Smith 240*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCSORSetSymmetric()`, `PCSORSetIterations()`, `PCEisenstatSetOmega()`, `PCSORSetOmega()` 241c60c7ad4SBarry Smith @*/ 242d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORGetOmega(PC pc, PetscReal *omega) 243d71ae5a4SJacob Faibussowitsch { 244c60c7ad4SBarry Smith PetscFunctionBegin; 245c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 246cac4c232SBarry Smith PetscUseMethod(pc, "PCSORGetOmega_C", (PC, PetscReal *), (pc, omega)); 2473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 248c60c7ad4SBarry Smith } 249c60c7ad4SBarry Smith 250c60c7ad4SBarry Smith /*@ 251c60c7ad4SBarry Smith PCSORGetIterations - Gets the number of inner iterations to 252c60c7ad4SBarry Smith be used by the SOR preconditioner. The default is 1. 253c60c7ad4SBarry Smith 254c3339decSBarry Smith Logically Collective 255c60c7ad4SBarry Smith 256c60c7ad4SBarry Smith Input Parameter: 257c60c7ad4SBarry Smith . pc - the preconditioner context 258c60c7ad4SBarry Smith 259d8d19677SJose E. Roman Output Parameters: 260c60c7ad4SBarry Smith + lits - number of local iterations, smoothings over just variables on processor 261c60c7ad4SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 262c60c7ad4SBarry Smith 263f1580f4eSBarry Smith Options Database Keys: 264c60c7ad4SBarry Smith + -pc_sor_its <its> - Sets number of iterations 265c60c7ad4SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 266c60c7ad4SBarry Smith 267c60c7ad4SBarry Smith Level: intermediate 268c60c7ad4SBarry Smith 269f1580f4eSBarry Smith Note: 27095452b02SPatrick Sanan When run on one processor the number of smoothings is lits*its 271c60c7ad4SBarry Smith 272*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCSORSetOmega()`, `PCSORSetSymmetric()`, `PCSORSetIterations()` 273c60c7ad4SBarry Smith @*/ 274d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORGetIterations(PC pc, PetscInt *its, PetscInt *lits) 275d71ae5a4SJacob Faibussowitsch { 276c60c7ad4SBarry Smith PetscFunctionBegin; 277c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 278cac4c232SBarry Smith PetscUseMethod(pc, "PCSORGetIterations_C", (PC, PetscInt *, PetscInt *), (pc, its, lits)); 2793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 280c60c7ad4SBarry Smith } 281c60c7ad4SBarry Smith 2824b9ad928SBarry Smith /*@ 2834b9ad928SBarry Smith PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR), 2844b9ad928SBarry Smith backward, or forward relaxation. The local variants perform SOR on 2854b9ad928SBarry Smith each processor. By default forward relaxation is used. 2864b9ad928SBarry Smith 287c3339decSBarry Smith Logically Collective 2884b9ad928SBarry Smith 2894b9ad928SBarry Smith Input Parameters: 2904b9ad928SBarry Smith + pc - the preconditioner context 2914b9ad928SBarry Smith - flag - one of the following 2924b9ad928SBarry Smith .vb 2934b9ad928SBarry Smith SOR_FORWARD_SWEEP 2944b9ad928SBarry Smith SOR_BACKWARD_SWEEP 2954b9ad928SBarry Smith SOR_SYMMETRIC_SWEEP 2964b9ad928SBarry Smith SOR_LOCAL_FORWARD_SWEEP 2974b9ad928SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 2984b9ad928SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 2994b9ad928SBarry Smith .ve 3004b9ad928SBarry Smith 3014b9ad928SBarry Smith Options Database Keys: 3024b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 3034b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 3044b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 3054b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 3064b9ad928SBarry Smith - -pc_sor_local_backward - Activates local backward version 3074b9ad928SBarry Smith 308f1580f4eSBarry Smith Note: 3094b9ad928SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 3104b9ad928SBarry Smith which can be chosen with the option 3114b9ad928SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 3124b9ad928SBarry Smith 3134b9ad928SBarry Smith Level: intermediate 3144b9ad928SBarry Smith 315*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCEisenstatSetOmega()`, `PCSORSetIterations()`, `PCSORSetOmega()` 3164b9ad928SBarry Smith @*/ 317d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORSetSymmetric(PC pc, MatSORType flag) 318d71ae5a4SJacob Faibussowitsch { 3194b9ad928SBarry Smith PetscFunctionBegin; 3200700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 321c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, flag, 2); 322cac4c232SBarry Smith PetscTryMethod(pc, "PCSORSetSymmetric_C", (PC, MatSORType), (pc, flag)); 3233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3244b9ad928SBarry Smith } 3254b9ad928SBarry Smith 3264b9ad928SBarry Smith /*@ 3274b9ad928SBarry Smith PCSORSetOmega - Sets the SOR relaxation coefficient, omega 3284b9ad928SBarry Smith (where omega = 1.0 by default). 3294b9ad928SBarry Smith 330c3339decSBarry Smith Logically Collective 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith Input Parameters: 3334b9ad928SBarry Smith + pc - the preconditioner context 3344b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2). 3354b9ad928SBarry Smith 3364b9ad928SBarry Smith Options Database Key: 3374b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith Level: intermediate 3404b9ad928SBarry Smith 341485168efSMatthew Knepley Note: 342f1580f4eSBarry Smith If omega != 1, you will need to set the `MAT_USE_INODE`S option to `PETSC_FALSE` on the matrix. 343485168efSMatthew Knepley 344*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCSORSetSymmetric()`, `PCSORSetIterations()`, `PCEisenstatSetOmega()`, `MatSetOption()` 3454b9ad928SBarry Smith @*/ 346d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORSetOmega(PC pc, PetscReal omega) 347d71ae5a4SJacob Faibussowitsch { 3484b9ad928SBarry Smith PetscFunctionBegin; 349c5eb9154SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 350c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc, omega, 2); 351cac4c232SBarry Smith PetscTryMethod(pc, "PCSORSetOmega_C", (PC, PetscReal), (pc, omega)); 3523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3534b9ad928SBarry Smith } 3544b9ad928SBarry Smith 3554b9ad928SBarry Smith /*@ 3564b9ad928SBarry Smith PCSORSetIterations - Sets the number of inner iterations to 3574b9ad928SBarry Smith be used by the SOR preconditioner. The default is 1. 3584b9ad928SBarry Smith 359c3339decSBarry Smith Logically Collective 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith Input Parameters: 3624b9ad928SBarry Smith + pc - the preconditioner context 3634b9ad928SBarry Smith . lits - number of local iterations, smoothings over just variables on processor 3644b9ad928SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 3654b9ad928SBarry Smith 366f1580f4eSBarry Smith Options Database Keys: 3674b9ad928SBarry Smith + -pc_sor_its <its> - Sets number of iterations 3684b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 3694b9ad928SBarry Smith 3704b9ad928SBarry Smith Level: intermediate 3714b9ad928SBarry Smith 372f1580f4eSBarry Smith Note: 37395452b02SPatrick Sanan When run on one processor the number of smoothings is lits*its 3744b9ad928SBarry Smith 375*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCSOR`, `PCSORSetOmega()`, `PCSORSetSymmetric()` 3764b9ad928SBarry Smith @*/ 377d71ae5a4SJacob Faibussowitsch PetscErrorCode PCSORSetIterations(PC pc, PetscInt its, PetscInt lits) 378d71ae5a4SJacob Faibussowitsch { 3794b9ad928SBarry Smith PetscFunctionBegin; 3800700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 381c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, its, 2); 382cac4c232SBarry Smith PetscTryMethod(pc, "PCSORSetIterations_C", (PC, PetscInt, PetscInt), (pc, its, lits)); 3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3844b9ad928SBarry Smith } 3854b9ad928SBarry Smith 3864b9ad928SBarry Smith /*MC 3874b9ad928SBarry Smith PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning 3884b9ad928SBarry Smith 3894b9ad928SBarry Smith Options Database Keys: 3904b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 3914b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 392a9510f2eSBarry Smith . -pc_sor_forward - Activates forward version 3934b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 394a9510f2eSBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version (default version) 3954b9ad928SBarry Smith . -pc_sor_local_backward - Activates local backward version 3964b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 397422a814eSBarry Smith . -pc_sor_diagonal_shift <shift> - shift the diagonal entries; useful if the matrix has zeros on the diagonal 398a9510f2eSBarry Smith . -pc_sor_its <its> - Sets number of iterations (default 1) 399a9510f2eSBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations (default 1) 4004b9ad928SBarry Smith 4014b9ad928SBarry Smith Level: beginner 4024b9ad928SBarry Smith 40395452b02SPatrick Sanan Notes: 404f1580f4eSBarry Smith Only implemented for the `MATAIJ` and `MATSEQBAIJ` matrix formats. 405f1580f4eSBarry Smith 4064b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 4074b9ad928SBarry Smith Jacobi with SOR on each block. 4084b9ad928SBarry Smith 409f1580f4eSBarry Smith For `MATAIJ` matrices if a diagonal entry is zero (and the diagonal shift is zero) then by default the inverse of that 41048773899SPierre Jolivet zero will be used and hence the `KSPSolve()` will terminate with `KSP_DIVERGED_NANORINF`. If the option 411f1580f4eSBarry Smith `KSPSetErrorIfNotConverged()` or -ksp_error_if_not_converged the code will terminate as soon as it detects the 412422a814eSBarry Smith zero pivot. 413422a814eSBarry Smith 414f1580f4eSBarry Smith For `MATSEQBAIJ` matrices this implements point-block SOR, but the omega, its, lits options are not supported. 41537a17b4dSBarry Smith 416f1580f4eSBarry Smith For `MATSEQBAIJ` the diagonal blocks are inverted using dense LU with partial pivoting. If a zero pivot is detected 417422a814eSBarry Smith the computation is stopped with an error 418422a814eSBarry Smith 419f1580f4eSBarry Smith If used with `KSPRICHARDSON` and no monitors the convergence test is skipped to improve speed, thus it always iterates 420f1580f4eSBarry Smith the maximum number of iterations you've selected for `KSP`. It is usually used in this mode as a smoother for multigrid. 42167991ba4SBarry Smith 422f1580f4eSBarry Smith If omega != 1, you will need to set the `MAT_USE_INODES` option to `PETSC_FALSE` on the matrix. 423485168efSMatthew Knepley 424*562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCJACOBI`, 425db781477SPatrick Sanan `PCSORSetIterations()`, `PCSORSetSymmetric()`, `PCSORSetOmega()`, `PCEISENSTAT`, `MatSetOption()` 4264b9ad928SBarry Smith M*/ 4274b9ad928SBarry Smith 428d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_SOR(PC pc) 429d71ae5a4SJacob Faibussowitsch { 4304b9ad928SBarry Smith PC_SOR *jac; 4314b9ad928SBarry Smith 4324b9ad928SBarry Smith PetscFunctionBegin; 4334dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&jac)); 4344b9ad928SBarry Smith 4354b9ad928SBarry Smith pc->ops->apply = PCApply_SOR; 4369d2471e0SBarry Smith pc->ops->applytranspose = PCApplyTranspose_SOR; 4374b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_SOR; 4384b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SOR; 4390a545947SLisandro Dalcin pc->ops->setup = NULL; 4404b9ad928SBarry Smith pc->ops->view = PCView_SOR; 4414b9ad928SBarry Smith pc->ops->destroy = PCDestroy_SOR; 4424b9ad928SBarry Smith pc->data = (void *)jac; 443d9bc8e36SBarry Smith jac->sym = SOR_LOCAL_SYMMETRIC_SWEEP; 4444b9ad928SBarry Smith jac->omega = 1.0; 44596fc60bcSBarry Smith jac->fshift = 0.0; 4464b9ad928SBarry Smith jac->its = 1; 4474b9ad928SBarry Smith jac->lits = 1; 4484b9ad928SBarry Smith 4499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetSymmetric_C", PCSORSetSymmetric_SOR)); 4509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetOmega_C", PCSORSetOmega_SOR)); 4519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORSetIterations_C", PCSORSetIterations_SOR)); 4529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetSymmetric_C", PCSORGetSymmetric_SOR)); 4539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetOmega_C", PCSORGetOmega_SOR)); 4549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCSORGetIterations_C", PCSORGetIterations_SOR)); 4553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4564b9ad928SBarry Smith } 457