1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines a Eisenstat trick SSOR preconditioner. This uses about 44b9ad928SBarry Smith %50 of the usual amount of floating point ops used for SSOR + Krylov 54b9ad928SBarry Smith method. But it requires actually solving the preconditioned problem 64b9ad928SBarry Smith with both left and right preconditioning. 74b9ad928SBarry Smith */ 8af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 94b9ad928SBarry Smith 104b9ad928SBarry Smith typedef struct { 114b9ad928SBarry Smith Mat shell,A; 1278c391d7SBarry Smith Vec b[2],diag; /* temporary storage for true right hand side */ 134b9ad928SBarry Smith PetscReal omega; 14ace3abfcSBarry Smith PetscBool usediag; /* indicates preconditioner should include diagonal scaling*/ 154b9ad928SBarry Smith } PC_Eisenstat; 164b9ad928SBarry Smith 176849ba73SBarry Smith static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x) 184b9ad928SBarry Smith { 194b9ad928SBarry Smith PC pc; 204b9ad928SBarry Smith PC_Eisenstat *eis; 214b9ad928SBarry Smith 224b9ad928SBarry Smith PetscFunctionBegin; 239566063dSJacob Faibussowitsch PetscCall(MatShellGetContext(mat,&pc)); 244b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 259566063dSJacob Faibussowitsch PetscCall(MatSOR(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x)); 264b9ad928SBarry Smith PetscFunctionReturn(0); 274b9ad928SBarry Smith } 284b9ad928SBarry Smith 296849ba73SBarry Smith static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y) 304b9ad928SBarry Smith { 314b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 32ace3abfcSBarry Smith PetscBool hasop; 334b9ad928SBarry Smith 344b9ad928SBarry Smith PetscFunctionBegin; 3589c6957cSBarry Smith if (eis->usediag) { 369566063dSJacob Faibussowitsch PetscCall(MatHasOperation(pc->pmat,MATOP_MULT_DIAGONAL_BLOCK,&hasop)); 3789c6957cSBarry Smith if (hasop) { 389566063dSJacob Faibussowitsch PetscCall(MatMultDiagonalBlock(pc->pmat,x,y)); 3989c6957cSBarry Smith } else { 409566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(y,x,eis->diag)); 4189c6957cSBarry Smith } 429566063dSJacob Faibussowitsch } else PetscCall(VecCopy(x,y)); 434b9ad928SBarry Smith PetscFunctionReturn(0); 444b9ad928SBarry Smith } 454b9ad928SBarry Smith 46946967b8SBarry Smith static PetscErrorCode PCPreSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x) 474b9ad928SBarry Smith { 484b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 49ace3abfcSBarry Smith PetscBool nonzero; 504b9ad928SBarry Smith 514b9ad928SBarry Smith PetscFunctionBegin; 5278c391d7SBarry Smith if (pc->presolvedone < 2) { 5308401ef6SPierre Jolivet PetscCheck(pc->mat == pc->pmat,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot have different mat and pmat"); 544b9ad928SBarry Smith /* swap shell matrix and true matrix */ 554b9ad928SBarry Smith eis->A = pc->mat; 564b9ad928SBarry Smith pc->mat = eis->shell; 574b9ad928SBarry Smith } 584b9ad928SBarry Smith 5978c391d7SBarry Smith if (!eis->b[pc->presolvedone-1]) { 609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b,&eis->b[pc->presolvedone-1])); 619566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->b[pc->presolvedone-1])); 6278c391d7SBarry Smith } 634b9ad928SBarry Smith 644b9ad928SBarry Smith /* if nonzero initial guess, modify x */ 659566063dSJacob Faibussowitsch PetscCall(KSPGetInitialGuessNonzero(ksp,&nonzero)); 664b9ad928SBarry Smith if (nonzero) { 679566063dSJacob Faibussowitsch PetscCall(VecCopy(x,eis->b[pc->presolvedone-1])); 689566063dSJacob Faibussowitsch PetscCall(MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,SOR_APPLY_UPPER,0.0,1,1,x)); 694b9ad928SBarry Smith } 704b9ad928SBarry Smith 71121471adSBarry Smith /* save true b, other option is to swap pointers */ 729566063dSJacob Faibussowitsch PetscCall(VecCopy(b,eis->b[pc->presolvedone-1])); 73121471adSBarry Smith 745c99c7daSBarry Smith /* modify b by (L + D/omega)^{-1} */ 759566063dSJacob Faibussowitsch PetscCall(MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),0.0,1,1,b)); 764b9ad928SBarry Smith PetscFunctionReturn(0); 774b9ad928SBarry Smith } 784b9ad928SBarry Smith 79946967b8SBarry Smith static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x) 804b9ad928SBarry Smith { 814b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 824b9ad928SBarry Smith 834b9ad928SBarry Smith PetscFunctionBegin; 844b9ad928SBarry Smith /* get back true b */ 859566063dSJacob Faibussowitsch PetscCall(VecCopy(eis->b[pc->presolvedone],b)); 86121471adSBarry Smith 87121471adSBarry Smith /* modify x by (U + D/omega)^{-1} */ 889566063dSJacob Faibussowitsch PetscCall(VecCopy(x,eis->b[pc->presolvedone])); 899566063dSJacob Faibussowitsch PetscCall(MatSOR(eis->A,eis->b[pc->presolvedone],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x)); 902fa5cd67SKarl Rupp if (!pc->presolvedone) pc->mat = eis->A; 914b9ad928SBarry Smith PetscFunctionReturn(0); 924b9ad928SBarry Smith } 934b9ad928SBarry Smith 9469d2c0f9SBarry Smith static PetscErrorCode PCReset_Eisenstat(PC pc) 954b9ad928SBarry Smith { 964b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 974b9ad928SBarry Smith 984b9ad928SBarry Smith PetscFunctionBegin; 999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eis->b[0])); 1009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eis->b[1])); 1019566063dSJacob Faibussowitsch PetscCall(MatDestroy(&eis->shell)); 1029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&eis->diag)); 10369d2c0f9SBarry Smith PetscFunctionReturn(0); 10469d2c0f9SBarry Smith } 10569d2c0f9SBarry Smith 10669d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Eisenstat(PC pc) 10769d2c0f9SBarry Smith { 10869d2c0f9SBarry Smith PetscFunctionBegin; 1099566063dSJacob Faibussowitsch PetscCall(PCReset_Eisenstat(pc)); 110*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetOmega_C",NULL)); 111*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetNoDiagonalScaling_C",NULL)); 112*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetOmega_C",NULL)); 113*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetNoDiagonalScaling_C",NULL)); 114*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",NULL)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 1164b9ad928SBarry Smith PetscFunctionReturn(0); 1174b9ad928SBarry Smith } 1184b9ad928SBarry Smith 1194416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Eisenstat(PetscOptionItems *PetscOptionsObject,PC pc) 1204b9ad928SBarry Smith { 1214b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 1228afaa268SBarry Smith PetscBool set,flg; 1234b9ad928SBarry Smith 1244b9ad928SBarry Smith PetscFunctionBegin; 125d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Eisenstat SSOR options"); 1269566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,NULL)); 1279566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatSetNoDiagonalScaling",eis->usediag ? PETSC_FALSE : PETSC_TRUE,&flg,&set)); 128c60c7ad4SBarry Smith if (set) { 1299566063dSJacob Faibussowitsch PetscCall(PCEisenstatSetNoDiagonalScaling(pc,flg)); 1304b9ad928SBarry Smith } 131d0609cedSBarry Smith PetscOptionsHeadEnd(); 1324b9ad928SBarry Smith PetscFunctionReturn(0); 1334b9ad928SBarry Smith } 1344b9ad928SBarry Smith 1356849ba73SBarry Smith static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer) 1364b9ad928SBarry Smith { 1374b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 138ace3abfcSBarry Smith PetscBool iascii; 1394b9ad928SBarry Smith 1404b9ad928SBarry Smith PetscFunctionBegin; 1419566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 14232077d6dSBarry Smith if (iascii) { 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," omega = %g\n",(double)eis->omega)); 1444b9ad928SBarry Smith if (eis->usediag) { 1459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using diagonal scaling (default)\n")); 1464b9ad928SBarry Smith } else { 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Not using diagonal scaling\n")); 1484b9ad928SBarry Smith } 1494b9ad928SBarry Smith } 1504b9ad928SBarry Smith PetscFunctionReturn(0); 1514b9ad928SBarry Smith } 1524b9ad928SBarry Smith 1536849ba73SBarry Smith static PetscErrorCode PCSetUp_Eisenstat(PC pc) 1544b9ad928SBarry Smith { 15513f74950SBarry Smith PetscInt M,N,m,n; 1564b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 1574b9ad928SBarry Smith 1584b9ad928SBarry Smith PetscFunctionBegin; 1594b9ad928SBarry Smith if (!pc->setupcalled) { 1609566063dSJacob Faibussowitsch PetscCall(MatGetSize(pc->mat,&M,&N)); 1619566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->mat,&m,&n)); 1629566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)pc),&eis->shell)); 1639566063dSJacob Faibussowitsch PetscCall(MatSetSizes(eis->shell,m,n,M,N)); 1649566063dSJacob Faibussowitsch PetscCall(MatSetType(eis->shell,MATSHELL)); 1659566063dSJacob Faibussowitsch PetscCall(MatSetUp(eis->shell)); 1669566063dSJacob Faibussowitsch PetscCall(MatShellSetContext(eis->shell,pc)); 1679566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->shell)); 1689566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(eis->shell,MATOP_MULT,(void (*)(void))PCMult_Eisenstat)); 1694b9ad928SBarry Smith } 1704b9ad928SBarry Smith if (!eis->usediag) PetscFunctionReturn(0); 1714b9ad928SBarry Smith if (!pc->setupcalled) { 1729566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat,&eis->diag,NULL)); 1739566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->diag)); 1744b9ad928SBarry Smith } 1759566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(pc->pmat,eis->diag)); 1764b9ad928SBarry Smith PetscFunctionReturn(0); 1774b9ad928SBarry Smith } 1784b9ad928SBarry Smith 1794b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 1804b9ad928SBarry Smith 1811e6b0712SBarry Smith static PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega) 1824b9ad928SBarry Smith { 183c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 1844b9ad928SBarry Smith 1854b9ad928SBarry Smith PetscFunctionBegin; 1862472a847SBarry Smith PetscCheck(omega > 0.0 && omega < 2.0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 1874b9ad928SBarry Smith eis->omega = omega; 1884b9ad928SBarry Smith PetscFunctionReturn(0); 1894b9ad928SBarry Smith } 1904b9ad928SBarry Smith 191c60c7ad4SBarry Smith static PetscErrorCode PCEisenstatSetNoDiagonalScaling_Eisenstat(PC pc,PetscBool flg) 1924b9ad928SBarry Smith { 193c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 1944b9ad928SBarry Smith 1954b9ad928SBarry Smith PetscFunctionBegin; 196c60c7ad4SBarry Smith eis->usediag = flg; 197c60c7ad4SBarry Smith PetscFunctionReturn(0); 198c60c7ad4SBarry Smith } 199c60c7ad4SBarry Smith 200c60c7ad4SBarry Smith static PetscErrorCode PCEisenstatGetOmega_Eisenstat(PC pc,PetscReal *omega) 201c60c7ad4SBarry Smith { 202c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 203c60c7ad4SBarry Smith 204c60c7ad4SBarry Smith PetscFunctionBegin; 205c60c7ad4SBarry Smith *omega = eis->omega; 206c60c7ad4SBarry Smith PetscFunctionReturn(0); 207c60c7ad4SBarry Smith } 208c60c7ad4SBarry Smith 209c60c7ad4SBarry Smith static PetscErrorCode PCEisenstatGetNoDiagonalScaling_Eisenstat(PC pc,PetscBool *flg) 210c60c7ad4SBarry Smith { 211c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 212c60c7ad4SBarry Smith 213c60c7ad4SBarry Smith PetscFunctionBegin; 214c60c7ad4SBarry Smith *flg = eis->usediag; 2154b9ad928SBarry Smith PetscFunctionReturn(0); 2164b9ad928SBarry Smith } 2174b9ad928SBarry Smith 2184b9ad928SBarry Smith /*@ 2194b9ad928SBarry Smith PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega, 2204b9ad928SBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default). 2214b9ad928SBarry Smith 222ad4df100SBarry Smith Logically Collective on PC 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith Input Parameters: 2254b9ad928SBarry Smith + pc - the preconditioner context 2264b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2) 2274b9ad928SBarry Smith 2284b9ad928SBarry Smith Options Database Key: 2294b9ad928SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 2304b9ad928SBarry Smith 2314b9ad928SBarry Smith Notes: 2324b9ad928SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 2334b9ad928SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 2344b9ad928SBarry Smith however, the preconditioned problem must be solved with both left 2354b9ad928SBarry Smith and right preconditioning. 2364b9ad928SBarry Smith 2374b9ad928SBarry Smith To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 2384b9ad928SBarry Smith which can be chosen with the database options 2394b9ad928SBarry Smith $ -pc_type sor -pc_sor_symmetric 2404b9ad928SBarry Smith 2414b9ad928SBarry Smith Level: intermediate 2424b9ad928SBarry Smith 243db781477SPatrick Sanan .seealso: `PCSORSetOmega()` 2444b9ad928SBarry Smith @*/ 2457087cfbeSBarry Smith PetscErrorCode PCEisenstatSetOmega(PC pc,PetscReal omega) 2464b9ad928SBarry Smith { 2474b9ad928SBarry Smith PetscFunctionBegin; 2480700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 249c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc,omega,2); 250cac4c232SBarry Smith PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega)); 2514b9ad928SBarry Smith PetscFunctionReturn(0); 2524b9ad928SBarry Smith } 2534b9ad928SBarry Smith 2544b9ad928SBarry Smith /*@ 255c60c7ad4SBarry Smith PCEisenstatSetNoDiagonalScaling - Causes the Eisenstat preconditioner 2564b9ad928SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 2574b9ad928SBarry Smith along the diagonal, this may save a small amount of work. 2584b9ad928SBarry Smith 259ad4df100SBarry Smith Logically Collective on PC 2604b9ad928SBarry Smith 261c60c7ad4SBarry Smith Input Parameters: 262c60c7ad4SBarry Smith + pc - the preconditioner context 263c60c7ad4SBarry Smith - flg - PETSC_TRUE turns off diagonal scaling inside the algorithm 2644b9ad928SBarry Smith 2654b9ad928SBarry Smith Options Database Key: 266c60c7ad4SBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling() 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith Level: intermediate 2694b9ad928SBarry Smith 2704b9ad928SBarry Smith Note: 271f9ff08acSPierre Jolivet If you use the KSPSetDiagonalScaling() or -ksp_diagonal_scale option then you will 2726aad120cSJose E. Roman likely want to use this routine since it will save you some unneeded flops. 2734b9ad928SBarry Smith 274db781477SPatrick Sanan .seealso: `PCEisenstatSetOmega()` 2754b9ad928SBarry Smith @*/ 276c60c7ad4SBarry Smith PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC pc,PetscBool flg) 2774b9ad928SBarry Smith { 2784b9ad928SBarry Smith PetscFunctionBegin; 2790700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 280cac4c232SBarry Smith PetscTryMethod(pc,"PCEisenstatSetNoDiagonalScaling_C",(PC,PetscBool),(pc,flg)); 281c60c7ad4SBarry Smith PetscFunctionReturn(0); 282c60c7ad4SBarry Smith } 283c60c7ad4SBarry Smith 284c60c7ad4SBarry Smith /*@ 285c60c7ad4SBarry Smith PCEisenstatGetOmega - Gets the SSOR relaxation coefficient, omega, 286c60c7ad4SBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default). 287c60c7ad4SBarry Smith 288c60c7ad4SBarry Smith Logically Collective on PC 289c60c7ad4SBarry Smith 290c60c7ad4SBarry Smith Input Parameter: 291c60c7ad4SBarry Smith . pc - the preconditioner context 292c60c7ad4SBarry Smith 293c60c7ad4SBarry Smith Output Parameter: 294c60c7ad4SBarry Smith . omega - relaxation coefficient (0 < omega < 2) 295c60c7ad4SBarry Smith 296c60c7ad4SBarry Smith Options Database Key: 297c60c7ad4SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 298c60c7ad4SBarry Smith 299c60c7ad4SBarry Smith Notes: 300c60c7ad4SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 301c60c7ad4SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 302c60c7ad4SBarry Smith however, the preconditioned problem must be solved with both left 303c60c7ad4SBarry Smith and right preconditioning. 304c60c7ad4SBarry Smith 305c60c7ad4SBarry Smith To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 306c60c7ad4SBarry Smith which can be chosen with the database options 307c60c7ad4SBarry Smith $ -pc_type sor -pc_sor_symmetric 308c60c7ad4SBarry Smith 309c60c7ad4SBarry Smith Level: intermediate 310c60c7ad4SBarry Smith 311db781477SPatrick Sanan .seealso: `PCSORGetOmega()`, `PCEisenstatSetOmega()` 312c60c7ad4SBarry Smith @*/ 313c60c7ad4SBarry Smith PetscErrorCode PCEisenstatGetOmega(PC pc,PetscReal *omega) 314c60c7ad4SBarry Smith { 315c60c7ad4SBarry Smith PetscFunctionBegin; 316c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 317cac4c232SBarry Smith PetscUseMethod(pc,"PCEisenstatGetOmega_C",(PC,PetscReal*),(pc,omega)); 318c60c7ad4SBarry Smith PetscFunctionReturn(0); 319c60c7ad4SBarry Smith } 320c60c7ad4SBarry Smith 321c60c7ad4SBarry Smith /*@ 322163d334eSBarry Smith PCEisenstatGetNoDiagonalScaling - Tells if the Eisenstat preconditioner 323c60c7ad4SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 324c60c7ad4SBarry Smith along the diagonal, this may save a small amount of work. 325c60c7ad4SBarry Smith 326c60c7ad4SBarry Smith Logically Collective on PC 327c60c7ad4SBarry Smith 328c60c7ad4SBarry Smith Input Parameter: 329c60c7ad4SBarry Smith . pc - the preconditioner context 330c60c7ad4SBarry Smith 331c60c7ad4SBarry Smith Output Parameter: 332c60c7ad4SBarry Smith . flg - PETSC_TRUE means there is no diagonal scaling applied 333c60c7ad4SBarry Smith 334c60c7ad4SBarry Smith Options Database Key: 335c60c7ad4SBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling() 336c60c7ad4SBarry Smith 337c60c7ad4SBarry Smith Level: intermediate 338c60c7ad4SBarry Smith 339c60c7ad4SBarry Smith Note: 340f9ff08acSPierre Jolivet If you use the KSPSetDiagonalScaling() or -ksp_diagonal_scale option then you will 3416aad120cSJose E. Roman likely want to use this routine since it will save you some unneeded flops. 342c60c7ad4SBarry Smith 343db781477SPatrick Sanan .seealso: `PCEisenstatGetOmega()` 344c60c7ad4SBarry Smith @*/ 345c60c7ad4SBarry Smith PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC pc,PetscBool *flg) 346c60c7ad4SBarry Smith { 347c60c7ad4SBarry Smith PetscFunctionBegin; 348c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 349cac4c232SBarry Smith PetscUseMethod(pc,"PCEisenstatGetNoDiagonalScaling_C",(PC,PetscBool*),(pc,flg)); 3504b9ad928SBarry Smith PetscFunctionReturn(0); 3514b9ad928SBarry Smith } 3524b9ad928SBarry Smith 3538066bbecSBarry Smith static PetscErrorCode PCPreSolveChangeRHS_Eisenstat(PC pc, PetscBool* change) 3548066bbecSBarry Smith { 3558066bbecSBarry Smith PetscFunctionBegin; 3568066bbecSBarry Smith *change = PETSC_TRUE; 3578066bbecSBarry Smith PetscFunctionReturn(0); 3588066bbecSBarry Smith } 3598066bbecSBarry Smith 3604b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 3614b9ad928SBarry Smith 3624b9ad928SBarry Smith /*MC 3634b9ad928SBarry Smith PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 3644b9ad928SBarry Smith preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 3654b9ad928SBarry Smith 3664b9ad928SBarry Smith Options Database Keys: 3674b9ad928SBarry Smith + -pc_eisenstat_omega <omega> - Sets omega 368c60c7ad4SBarry Smith - -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling() 3694b9ad928SBarry Smith 3704b9ad928SBarry Smith Level: beginner 3714b9ad928SBarry Smith 37295452b02SPatrick Sanan Notes: 37395452b02SPatrick Sanan Only implemented for the SeqAIJ matrix format. 3744b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 3754b9ad928SBarry Smith Jacobi with SOR on each block. 3764b9ad928SBarry Smith 377db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 378db781477SPatrick Sanan `PCEisenstatSetNoDiagonalScaling()`, `PCEisenstatSetOmega()`, `PCSOR` 3794b9ad928SBarry Smith M*/ 3804b9ad928SBarry Smith 3818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc) 3824b9ad928SBarry Smith { 3834b9ad928SBarry Smith PC_Eisenstat *eis; 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith PetscFunctionBegin; 3869566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&eis)); 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith pc->ops->apply = PCApply_Eisenstat; 389dc231df0SBarry Smith pc->ops->presolve = PCPreSolve_Eisenstat; 390dc231df0SBarry Smith pc->ops->postsolve = PCPostSolve_Eisenstat; 3910a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 3924b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 3934b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Eisenstat; 39469d2c0f9SBarry Smith pc->ops->reset = PCReset_Eisenstat; 3954b9ad928SBarry Smith pc->ops->view = PCView_Eisenstat; 3964b9ad928SBarry Smith pc->ops->setup = PCSetUp_Eisenstat; 3974b9ad928SBarry Smith 3983ec1f749SStefano Zampini pc->data = eis; 3994b9ad928SBarry Smith eis->omega = 1.0; 4000a545947SLisandro Dalcin eis->b[0] = NULL; 4010a545947SLisandro Dalcin eis->b[1] = NULL; 4020a545947SLisandro Dalcin eis->diag = NULL; 4034b9ad928SBarry Smith eis->usediag = PETSC_TRUE; 4044b9ad928SBarry Smith 4059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetOmega_C",PCEisenstatSetOmega_Eisenstat)); 4069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetNoDiagonalScaling_C",PCEisenstatSetNoDiagonalScaling_Eisenstat)); 4079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetOmega_C",PCEisenstatGetOmega_Eisenstat)); 4089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetNoDiagonalScaling_C",PCEisenstatGetNoDiagonalScaling_Eisenstat)); 4099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_Eisenstat)); 4104b9ad928SBarry Smith PetscFunctionReturn(0); 4114b9ad928SBarry Smith } 412