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 174b9ad928SBarry Smith 186849ba73SBarry Smith static PetscErrorCode PCMult_Eisenstat(Mat mat,Vec b,Vec x) 194b9ad928SBarry Smith { 20dfbe8321SBarry Smith PetscErrorCode ierr; 214b9ad928SBarry Smith PC pc; 224b9ad928SBarry Smith PC_Eisenstat *eis; 234b9ad928SBarry Smith 244b9ad928SBarry Smith PetscFunctionBegin; 254b9ad928SBarry Smith ierr = MatShellGetContext(mat,(void**)&pc);CHKERRQ(ierr); 264b9ad928SBarry Smith eis = (PC_Eisenstat*)pc->data; 2741f059aeSBarry Smith ierr = MatSOR(eis->A,b,eis->omega,SOR_EISENSTAT,0.0,1,1,x);CHKERRQ(ierr); 284b9ad928SBarry Smith PetscFunctionReturn(0); 294b9ad928SBarry Smith } 304b9ad928SBarry Smith 316849ba73SBarry Smith static PetscErrorCode PCApply_Eisenstat(PC pc,Vec x,Vec y) 324b9ad928SBarry Smith { 334b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 34dfbe8321SBarry Smith PetscErrorCode ierr; 35ace3abfcSBarry Smith PetscBool hasop; 364b9ad928SBarry Smith 374b9ad928SBarry Smith PetscFunctionBegin; 3889c6957cSBarry Smith if (eis->usediag) { 3989c6957cSBarry Smith ierr = MatHasOperation(pc->pmat,MATOP_MULT_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr); 4089c6957cSBarry Smith if (hasop) { 4189c6957cSBarry Smith ierr = MatMultDiagonalBlock(pc->pmat,x,y);CHKERRQ(ierr); 4289c6957cSBarry Smith } else { 4389c6957cSBarry Smith ierr = VecPointwiseMult(y,x,eis->diag);CHKERRQ(ierr); 4489c6957cSBarry Smith } 4589c6957cSBarry Smith } else {ierr = VecCopy(x,y);CHKERRQ(ierr);} 464b9ad928SBarry Smith PetscFunctionReturn(0); 474b9ad928SBarry Smith } 484b9ad928SBarry Smith 49946967b8SBarry Smith static PetscErrorCode PCPreSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x) 504b9ad928SBarry Smith { 514b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 52ace3abfcSBarry Smith PetscBool nonzero; 53dfbe8321SBarry Smith PetscErrorCode ierr; 544b9ad928SBarry Smith 554b9ad928SBarry Smith PetscFunctionBegin; 5678c391d7SBarry Smith if (pc->presolvedone < 2) { 57ce94432eSBarry Smith if (pc->mat != pc->pmat) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot have different mat and pmat"); 584b9ad928SBarry Smith /* swap shell matrix and true matrix */ 594b9ad928SBarry Smith eis->A = pc->mat; 604b9ad928SBarry Smith pc->mat = eis->shell; 614b9ad928SBarry Smith } 624b9ad928SBarry Smith 6378c391d7SBarry Smith if (!eis->b[pc->presolvedone-1]) { 6478c391d7SBarry Smith ierr = VecDuplicate(b,&eis->b[pc->presolvedone-1]);CHKERRQ(ierr); 653bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->b[pc->presolvedone-1]);CHKERRQ(ierr); 6678c391d7SBarry Smith } 674b9ad928SBarry Smith 684b9ad928SBarry Smith /* if nonzero initial guess, modify x */ 694b9ad928SBarry Smith ierr = KSPGetInitialGuessNonzero(ksp,&nonzero);CHKERRQ(ierr); 704b9ad928SBarry Smith if (nonzero) { 7178c391d7SBarry Smith ierr = VecCopy(x,eis->b[pc->presolvedone-1]);CHKERRQ(ierr); 7278c391d7SBarry Smith ierr = MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,SOR_APPLY_UPPER,0.0,1,1,x);CHKERRQ(ierr); 734b9ad928SBarry Smith } 744b9ad928SBarry Smith 75121471adSBarry Smith /* save true b, other option is to swap pointers */ 7678c391d7SBarry Smith ierr = VecCopy(b,eis->b[pc->presolvedone-1]);CHKERRQ(ierr); 77121471adSBarry Smith 785c99c7daSBarry Smith /* modify b by (L + D/omega)^{-1} */ 7978c391d7SBarry Smith ierr = MatSOR(eis->A,eis->b[pc->presolvedone-1],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_FORWARD_SWEEP),0.0,1,1,b);CHKERRQ(ierr); 804b9ad928SBarry Smith PetscFunctionReturn(0); 814b9ad928SBarry Smith } 824b9ad928SBarry Smith 83946967b8SBarry Smith static PetscErrorCode PCPostSolve_Eisenstat(PC pc,KSP ksp,Vec b,Vec x) 844b9ad928SBarry Smith { 854b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 86dfbe8321SBarry Smith PetscErrorCode ierr; 874b9ad928SBarry Smith 884b9ad928SBarry Smith PetscFunctionBegin; 894b9ad928SBarry Smith /* get back true b */ 9078c391d7SBarry Smith ierr = VecCopy(eis->b[pc->presolvedone],b);CHKERRQ(ierr); 91121471adSBarry Smith 92121471adSBarry Smith /* modify x by (U + D/omega)^{-1} */ 9378c391d7SBarry Smith ierr = VecCopy(x,eis->b[pc->presolvedone]);CHKERRQ(ierr); 9478c391d7SBarry Smith ierr = MatSOR(eis->A,eis->b[pc->presolvedone],eis->omega,(MatSORType)(SOR_ZERO_INITIAL_GUESS | SOR_LOCAL_BACKWARD_SWEEP),0.0,1,1,x);CHKERRQ(ierr); 952fa5cd67SKarl Rupp if (!pc->presolvedone) pc->mat = eis->A; 964b9ad928SBarry Smith PetscFunctionReturn(0); 974b9ad928SBarry Smith } 984b9ad928SBarry Smith 9969d2c0f9SBarry Smith static PetscErrorCode PCReset_Eisenstat(PC pc) 1004b9ad928SBarry Smith { 1014b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 102dfbe8321SBarry Smith PetscErrorCode ierr; 1034b9ad928SBarry Smith 1044b9ad928SBarry Smith PetscFunctionBegin; 10578c391d7SBarry Smith ierr = VecDestroy(&eis->b[0]);CHKERRQ(ierr); 10678c391d7SBarry Smith ierr = VecDestroy(&eis->b[1]);CHKERRQ(ierr); 1076bf464f9SBarry Smith ierr = MatDestroy(&eis->shell);CHKERRQ(ierr); 1086bf464f9SBarry Smith ierr = VecDestroy(&eis->diag);CHKERRQ(ierr); 10969d2c0f9SBarry Smith PetscFunctionReturn(0); 11069d2c0f9SBarry Smith } 11169d2c0f9SBarry Smith 11269d2c0f9SBarry Smith static PetscErrorCode PCDestroy_Eisenstat(PC pc) 11369d2c0f9SBarry Smith { 11469d2c0f9SBarry Smith PetscErrorCode ierr; 11569d2c0f9SBarry Smith 11669d2c0f9SBarry Smith PetscFunctionBegin; 1173f34a243SBarry Smith ierr = PCReset_Eisenstat(pc);CHKERRQ(ierr); 118c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 1194b9ad928SBarry Smith PetscFunctionReturn(0); 1204b9ad928SBarry Smith } 1214b9ad928SBarry Smith 1224416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Eisenstat(PetscOptionItems *PetscOptionsObject,PC pc) 1234b9ad928SBarry Smith { 1244b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 125dfbe8321SBarry Smith PetscErrorCode ierr; 1268afaa268SBarry Smith PetscBool set,flg; 1274b9ad928SBarry Smith 1284b9ad928SBarry Smith PetscFunctionBegin; 129e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Eisenstat SSOR options");CHKERRQ(ierr); 1300298fd71SBarry Smith ierr = PetscOptionsReal("-pc_eisenstat_omega","Relaxation factor 0 < omega < 2","PCEisenstatSetOmega",eis->omega,&eis->omega,NULL);CHKERRQ(ierr); 131c60c7ad4SBarry Smith ierr = PetscOptionsBool("-pc_eisenstat_no_diagonal_scaling","Do not use standard diagonal scaling","PCEisenstatSetNoDiagonalScaling",eis->usediag ? PETSC_FALSE : PETSC_TRUE,&flg,&set);CHKERRQ(ierr); 132c60c7ad4SBarry Smith if (set) { 133c60c7ad4SBarry Smith ierr = PCEisenstatSetNoDiagonalScaling(pc,flg);CHKERRQ(ierr); 1344b9ad928SBarry Smith } 1354b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 1364b9ad928SBarry Smith PetscFunctionReturn(0); 1374b9ad928SBarry Smith } 1384b9ad928SBarry Smith 1396849ba73SBarry Smith static PetscErrorCode PCView_Eisenstat(PC pc,PetscViewer viewer) 1404b9ad928SBarry Smith { 1414b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 142dfbe8321SBarry Smith PetscErrorCode ierr; 143ace3abfcSBarry Smith PetscBool iascii; 1444b9ad928SBarry Smith 1454b9ad928SBarry Smith PetscFunctionBegin; 146251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 14732077d6dSBarry Smith if (iascii) { 14857622a8eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: omega = %g\n",(double)eis->omega);CHKERRQ(ierr); 1494b9ad928SBarry Smith if (eis->usediag) { 1504b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Using diagonal scaling (default)\n");CHKERRQ(ierr); 1514b9ad928SBarry Smith } else { 1524b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Eisenstat: Not using diagonal scaling\n");CHKERRQ(ierr); 1534b9ad928SBarry Smith } 1544b9ad928SBarry Smith } 1554b9ad928SBarry Smith PetscFunctionReturn(0); 1564b9ad928SBarry Smith } 1574b9ad928SBarry Smith 1586849ba73SBarry Smith static PetscErrorCode PCSetUp_Eisenstat(PC pc) 1594b9ad928SBarry Smith { 160dfbe8321SBarry Smith PetscErrorCode ierr; 16113f74950SBarry Smith PetscInt M,N,m,n; 1624b9ad928SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 1634b9ad928SBarry Smith 1644b9ad928SBarry Smith PetscFunctionBegin; 1654b9ad928SBarry Smith if (!pc->setupcalled) { 1664b9ad928SBarry Smith ierr = MatGetSize(pc->mat,&M,&N);CHKERRQ(ierr); 1674b9ad928SBarry Smith ierr = MatGetLocalSize(pc->mat,&m,&n);CHKERRQ(ierr); 168ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)pc),&eis->shell);CHKERRQ(ierr); 16964aae45aSBarry Smith ierr = MatSetSizes(eis->shell,m,n,M,N);CHKERRQ(ierr); 170f204ca49SKris Buschelman ierr = MatSetType(eis->shell,MATSHELL);CHKERRQ(ierr); 1717ae8954aSSatish Balay ierr = MatSetUp(eis->shell);CHKERRQ(ierr); 172f204ca49SKris Buschelman ierr = MatShellSetContext(eis->shell,(void*)pc);CHKERRQ(ierr); 1733bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->shell);CHKERRQ(ierr); 1744b9ad928SBarry Smith ierr = MatShellSetOperation(eis->shell,MATOP_MULT,(void (*)(void))PCMult_Eisenstat);CHKERRQ(ierr); 1754b9ad928SBarry Smith } 1764b9ad928SBarry Smith if (!eis->usediag) PetscFunctionReturn(0); 1774b9ad928SBarry Smith if (!pc->setupcalled) { 1782a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&eis->diag,0);CHKERRQ(ierr); 1793bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)eis->diag);CHKERRQ(ierr); 1804b9ad928SBarry Smith } 1814b9ad928SBarry Smith ierr = MatGetDiagonal(pc->pmat,eis->diag);CHKERRQ(ierr); 1824b9ad928SBarry Smith PetscFunctionReturn(0); 1834b9ad928SBarry Smith } 1844b9ad928SBarry Smith 1854b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 1864b9ad928SBarry Smith 1871e6b0712SBarry Smith static PetscErrorCode PCEisenstatSetOmega_Eisenstat(PC pc,PetscReal omega) 1884b9ad928SBarry Smith { 189c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 1904b9ad928SBarry Smith 1914b9ad928SBarry Smith PetscFunctionBegin; 192ce94432eSBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 1934b9ad928SBarry Smith eis->omega = omega; 1944b9ad928SBarry Smith PetscFunctionReturn(0); 1954b9ad928SBarry Smith } 1964b9ad928SBarry Smith 197c60c7ad4SBarry Smith static PetscErrorCode PCEisenstatSetNoDiagonalScaling_Eisenstat(PC pc,PetscBool flg) 1984b9ad928SBarry Smith { 199c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 2004b9ad928SBarry Smith 2014b9ad928SBarry Smith PetscFunctionBegin; 202c60c7ad4SBarry Smith eis->usediag = flg; 203c60c7ad4SBarry Smith PetscFunctionReturn(0); 204c60c7ad4SBarry Smith } 205c60c7ad4SBarry Smith 206c60c7ad4SBarry Smith static PetscErrorCode PCEisenstatGetOmega_Eisenstat(PC pc,PetscReal *omega) 207c60c7ad4SBarry Smith { 208c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 209c60c7ad4SBarry Smith 210c60c7ad4SBarry Smith PetscFunctionBegin; 211c60c7ad4SBarry Smith *omega = eis->omega; 212c60c7ad4SBarry Smith PetscFunctionReturn(0); 213c60c7ad4SBarry Smith } 214c60c7ad4SBarry Smith 215c60c7ad4SBarry Smith static PetscErrorCode PCEisenstatGetNoDiagonalScaling_Eisenstat(PC pc,PetscBool *flg) 216c60c7ad4SBarry Smith { 217c60c7ad4SBarry Smith PC_Eisenstat *eis = (PC_Eisenstat*)pc->data; 218c60c7ad4SBarry Smith 219c60c7ad4SBarry Smith PetscFunctionBegin; 220c60c7ad4SBarry Smith *flg = eis->usediag; 2214b9ad928SBarry Smith PetscFunctionReturn(0); 2224b9ad928SBarry Smith } 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith /*@ 2254b9ad928SBarry Smith PCEisenstatSetOmega - Sets the SSOR relaxation coefficient, omega, 2264b9ad928SBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default). 2274b9ad928SBarry Smith 228ad4df100SBarry Smith Logically Collective on PC 2294b9ad928SBarry Smith 2304b9ad928SBarry Smith Input Parameters: 2314b9ad928SBarry Smith + pc - the preconditioner context 2324b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2) 2334b9ad928SBarry Smith 2344b9ad928SBarry Smith Options Database Key: 2354b9ad928SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 2364b9ad928SBarry Smith 2374b9ad928SBarry Smith Notes: 2384b9ad928SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 2394b9ad928SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 2404b9ad928SBarry Smith however, the preconditioned problem must be solved with both left 2414b9ad928SBarry Smith and right preconditioning. 2424b9ad928SBarry Smith 2434b9ad928SBarry Smith To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 2444b9ad928SBarry Smith which can be chosen with the database options 2454b9ad928SBarry Smith $ -pc_type sor -pc_sor_symmetric 2464b9ad928SBarry Smith 2474b9ad928SBarry Smith Level: intermediate 2484b9ad928SBarry Smith 2494b9ad928SBarry Smith .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega 2504b9ad928SBarry Smith 2514b9ad928SBarry Smith .seealso: PCSORSetOmega() 2524b9ad928SBarry Smith @*/ 2537087cfbeSBarry Smith PetscErrorCode PCEisenstatSetOmega(PC pc,PetscReal omega) 2544b9ad928SBarry Smith { 2554ac538c5SBarry Smith PetscErrorCode ierr; 2564b9ad928SBarry Smith 2574b9ad928SBarry Smith PetscFunctionBegin; 2580700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 259c5eb9154SBarry Smith PetscValidLogicalCollectiveReal(pc,omega,2); 2604ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCEisenstatSetOmega_C",(PC,PetscReal),(pc,omega));CHKERRQ(ierr); 2614b9ad928SBarry Smith PetscFunctionReturn(0); 2624b9ad928SBarry Smith } 2634b9ad928SBarry Smith 2644b9ad928SBarry Smith /*@ 265c60c7ad4SBarry Smith PCEisenstatSetNoDiagonalScaling - Causes the Eisenstat preconditioner 2664b9ad928SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 2674b9ad928SBarry Smith along the diagonal, this may save a small amount of work. 2684b9ad928SBarry Smith 269ad4df100SBarry Smith Logically Collective on PC 2704b9ad928SBarry Smith 271c60c7ad4SBarry Smith Input Parameters: 272c60c7ad4SBarry Smith + pc - the preconditioner context 273c60c7ad4SBarry Smith - flg - PETSC_TRUE turns off diagonal scaling inside the algorithm 2744b9ad928SBarry Smith 2754b9ad928SBarry Smith Options Database Key: 276c60c7ad4SBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling() 2774b9ad928SBarry Smith 2784b9ad928SBarry Smith Level: intermediate 2794b9ad928SBarry Smith 2804b9ad928SBarry Smith Note: 2814b9ad928SBarry Smith If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will 2824b9ad928SBarry Smith likley want to use this routine since it will save you some unneeded flops. 2834b9ad928SBarry Smith 2844b9ad928SBarry Smith .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR 2854b9ad928SBarry Smith 2864b9ad928SBarry Smith .seealso: PCEisenstatSetOmega() 2874b9ad928SBarry Smith @*/ 288c60c7ad4SBarry Smith PetscErrorCode PCEisenstatSetNoDiagonalScaling(PC pc,PetscBool flg) 2894b9ad928SBarry Smith { 2904ac538c5SBarry Smith PetscErrorCode ierr; 2914b9ad928SBarry Smith 2924b9ad928SBarry Smith PetscFunctionBegin; 2930700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 294c60c7ad4SBarry Smith ierr = PetscTryMethod(pc,"PCEisenstatSetNoDiagonalScaling_C",(PC,PetscBool),(pc,flg));CHKERRQ(ierr); 295c60c7ad4SBarry Smith PetscFunctionReturn(0); 296c60c7ad4SBarry Smith } 297c60c7ad4SBarry Smith 298c60c7ad4SBarry Smith /*@ 299c60c7ad4SBarry Smith PCEisenstatGetOmega - Gets the SSOR relaxation coefficient, omega, 300c60c7ad4SBarry Smith to use with Eisenstat's trick (where omega = 1.0 by default). 301c60c7ad4SBarry Smith 302c60c7ad4SBarry Smith Logically Collective on PC 303c60c7ad4SBarry Smith 304c60c7ad4SBarry Smith Input Parameter: 305c60c7ad4SBarry Smith . pc - the preconditioner context 306c60c7ad4SBarry Smith 307c60c7ad4SBarry Smith Output Parameter: 308c60c7ad4SBarry Smith . omega - relaxation coefficient (0 < omega < 2) 309c60c7ad4SBarry Smith 310c60c7ad4SBarry Smith Options Database Key: 311c60c7ad4SBarry Smith . -pc_eisenstat_omega <omega> - Sets omega 312c60c7ad4SBarry Smith 313c60c7ad4SBarry Smith Notes: 314c60c7ad4SBarry Smith The Eisenstat trick implementation of SSOR requires about 50% of the 315c60c7ad4SBarry Smith usual amount of floating point operations used for SSOR + Krylov method; 316c60c7ad4SBarry Smith however, the preconditioned problem must be solved with both left 317c60c7ad4SBarry Smith and right preconditioning. 318c60c7ad4SBarry Smith 319c60c7ad4SBarry Smith To use SSOR without the Eisenstat trick, employ the PCSOR preconditioner, 320c60c7ad4SBarry Smith which can be chosen with the database options 321c60c7ad4SBarry Smith $ -pc_type sor -pc_sor_symmetric 322c60c7ad4SBarry Smith 323c60c7ad4SBarry Smith Level: intermediate 324c60c7ad4SBarry Smith 325c60c7ad4SBarry Smith .keywords: PC, Eisenstat, set, SOR, SSOR, relaxation, omega 326c60c7ad4SBarry Smith 327c60c7ad4SBarry Smith .seealso: PCSORGetOmega(), PCEisenstatSetOmega() 328c60c7ad4SBarry Smith @*/ 329c60c7ad4SBarry Smith PetscErrorCode PCEisenstatGetOmega(PC pc,PetscReal *omega) 330c60c7ad4SBarry Smith { 331c60c7ad4SBarry Smith PetscErrorCode ierr; 332c60c7ad4SBarry Smith 333c60c7ad4SBarry Smith PetscFunctionBegin; 334c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 335c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCEisenstatGetOmega_C",(PC,PetscReal*),(pc,omega));CHKERRQ(ierr); 336c60c7ad4SBarry Smith PetscFunctionReturn(0); 337c60c7ad4SBarry Smith } 338c60c7ad4SBarry Smith 339c60c7ad4SBarry Smith /*@ 340163d334eSBarry Smith PCEisenstatGetNoDiagonalScaling - Tells if the Eisenstat preconditioner 341c60c7ad4SBarry Smith not to do additional diagonal preconditioning. For matrices with a constant 342c60c7ad4SBarry Smith along the diagonal, this may save a small amount of work. 343c60c7ad4SBarry Smith 344c60c7ad4SBarry Smith Logically Collective on PC 345c60c7ad4SBarry Smith 346c60c7ad4SBarry Smith Input Parameter: 347c60c7ad4SBarry Smith . pc - the preconditioner context 348c60c7ad4SBarry Smith 349c60c7ad4SBarry Smith Output Parameter: 350c60c7ad4SBarry Smith . flg - PETSC_TRUE means there is no diagonal scaling applied 351c60c7ad4SBarry Smith 352c60c7ad4SBarry Smith Options Database Key: 353c60c7ad4SBarry Smith . -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling() 354c60c7ad4SBarry Smith 355c60c7ad4SBarry Smith Level: intermediate 356c60c7ad4SBarry Smith 357c60c7ad4SBarry Smith Note: 358c60c7ad4SBarry Smith If you use the KPSSetDiagonalScaling() or -ksp_diagonal_scale option then you will 359c60c7ad4SBarry Smith likley want to use this routine since it will save you some unneeded flops. 360c60c7ad4SBarry Smith 361c60c7ad4SBarry Smith .keywords: PC, Eisenstat, use, diagonal, scaling, SSOR 362c60c7ad4SBarry Smith 363c60c7ad4SBarry Smith .seealso: PCEisenstatGetOmega() 364c60c7ad4SBarry Smith @*/ 365c60c7ad4SBarry Smith PetscErrorCode PCEisenstatGetNoDiagonalScaling(PC pc,PetscBool *flg) 366c60c7ad4SBarry Smith { 367c60c7ad4SBarry Smith PetscErrorCode ierr; 368c60c7ad4SBarry Smith 369c60c7ad4SBarry Smith PetscFunctionBegin; 370c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 371163d334eSBarry Smith ierr = PetscUseMethod(pc,"PCEisenstatGetNoDiagonalScaling_C",(PC,PetscBool*),(pc,flg));CHKERRQ(ierr); 3724b9ad928SBarry Smith PetscFunctionReturn(0); 3734b9ad928SBarry Smith } 3744b9ad928SBarry Smith 375*8066bbecSBarry Smith static PetscErrorCode PCPreSolveChangeRHS_Eisenstat(PC pc, PetscBool* change) 376*8066bbecSBarry Smith { 377*8066bbecSBarry Smith PetscFunctionBegin; 378*8066bbecSBarry Smith *change = PETSC_TRUE; 379*8066bbecSBarry Smith PetscFunctionReturn(0); 380*8066bbecSBarry Smith } 381*8066bbecSBarry Smith 3824b9ad928SBarry Smith /* --------------------------------------------------------------------*/ 3834b9ad928SBarry Smith 3844b9ad928SBarry Smith /*MC 3854b9ad928SBarry Smith PCEISENSTAT - An implementation of SSOR (symmetric successive over relaxation, symmetric Gauss-Seidel) 3864b9ad928SBarry Smith preconditioning that incorporates Eisenstat's trick to reduce the amount of computation needed. 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith Options Database Keys: 3894b9ad928SBarry Smith + -pc_eisenstat_omega <omega> - Sets omega 390c60c7ad4SBarry Smith - -pc_eisenstat_no_diagonal_scaling - Activates PCEisenstatSetNoDiagonalScaling() 3914b9ad928SBarry Smith 3924b9ad928SBarry Smith Level: beginner 3934b9ad928SBarry Smith 3944b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel, Eisenstat's trick 3954b9ad928SBarry Smith 3964b9ad928SBarry Smith Notes: Only implemented for the SeqAIJ matrix format. 3974b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 3984b9ad928SBarry Smith Jacobi with SOR on each block. 3994b9ad928SBarry Smith 4004b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 401c60c7ad4SBarry Smith PCEisenstatSetNoDiagonalScaling(), PCEisenstatSetOmega(), PCSOR 4024b9ad928SBarry Smith M*/ 4034b9ad928SBarry Smith 4048cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Eisenstat(PC pc) 4054b9ad928SBarry Smith { 406dfbe8321SBarry Smith PetscErrorCode ierr; 4074b9ad928SBarry Smith PC_Eisenstat *eis; 4084b9ad928SBarry Smith 4094b9ad928SBarry Smith PetscFunctionBegin; 410b00a9115SJed Brown ierr = PetscNewLog(pc,&eis);CHKERRQ(ierr); 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith pc->ops->apply = PCApply_Eisenstat; 413dc231df0SBarry Smith pc->ops->presolve = PCPreSolve_Eisenstat; 414dc231df0SBarry Smith pc->ops->postsolve = PCPostSolve_Eisenstat; 4154b9ad928SBarry Smith pc->ops->applyrichardson = 0; 4164b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Eisenstat; 4174b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Eisenstat; 41869d2c0f9SBarry Smith pc->ops->reset = PCReset_Eisenstat; 4194b9ad928SBarry Smith pc->ops->view = PCView_Eisenstat; 4204b9ad928SBarry Smith pc->ops->setup = PCSetUp_Eisenstat; 4214b9ad928SBarry Smith 4224b9ad928SBarry Smith pc->data = (void*)eis; 4234b9ad928SBarry Smith eis->omega = 1.0; 42478c391d7SBarry Smith eis->b[0] = 0; 42578c391d7SBarry Smith eis->b[1] = 0; 4264b9ad928SBarry Smith eis->diag = 0; 4274b9ad928SBarry Smith eis->usediag = PETSC_TRUE; 4284b9ad928SBarry Smith 429bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetOmega_C",PCEisenstatSetOmega_Eisenstat);CHKERRQ(ierr); 430c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatSetNoDiagonalScaling_C",PCEisenstatSetNoDiagonalScaling_Eisenstat);CHKERRQ(ierr); 431c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetOmega_C",PCEisenstatGetOmega_Eisenstat);CHKERRQ(ierr); 432c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCEisenstatGetNoDiagonalScaling_C",PCEisenstatGetNoDiagonalScaling_Eisenstat);CHKERRQ(ierr); 433*8066bbecSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCPreSolveChangeRHS_C",PCPreSolveChangeRHS_Eisenstat);CHKERRQ(ierr); 4344b9ad928SBarry Smith PetscFunctionReturn(0); 4354b9ad928SBarry Smith } 436