1*4b9ad928SBarry Smith /*$Id: sor.c,v 1.103 2001/08/21 21:03:17 bsmith Exp $*/ 2*4b9ad928SBarry Smith 3*4b9ad928SBarry Smith /* 4*4b9ad928SBarry Smith Defines a (S)SOR preconditioner for any Mat implementation 5*4b9ad928SBarry Smith */ 6*4b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 7*4b9ad928SBarry Smith 8*4b9ad928SBarry Smith typedef struct { 9*4b9ad928SBarry Smith int its; /* inner iterations, number of sweeps */ 10*4b9ad928SBarry Smith int lits; /* local inner iterations, number of sweeps applied by the local matrix mat->A */ 11*4b9ad928SBarry Smith MatSORType sym; /* forward, reverse, symmetric etc. */ 12*4b9ad928SBarry Smith PetscReal omega; 13*4b9ad928SBarry Smith } PC_SOR; 14*4b9ad928SBarry Smith 15*4b9ad928SBarry Smith #undef __FUNCT__ 16*4b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_SOR" 17*4b9ad928SBarry Smith static int PCDestroy_SOR(PC pc) 18*4b9ad928SBarry Smith { 19*4b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 20*4b9ad928SBarry Smith int ierr; 21*4b9ad928SBarry Smith 22*4b9ad928SBarry Smith PetscFunctionBegin; 23*4b9ad928SBarry Smith ierr = PetscFree(jac);CHKERRQ(ierr); 24*4b9ad928SBarry Smith PetscFunctionReturn(0); 25*4b9ad928SBarry Smith } 26*4b9ad928SBarry Smith 27*4b9ad928SBarry Smith #undef __FUNCT__ 28*4b9ad928SBarry Smith #define __FUNCT__ "PCApply_SOR" 29*4b9ad928SBarry Smith static int PCApply_SOR(PC pc,Vec x,Vec y) 30*4b9ad928SBarry Smith { 31*4b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 32*4b9ad928SBarry Smith int ierr,flag = jac->sym | SOR_ZERO_INITIAL_GUESS; 33*4b9ad928SBarry Smith 34*4b9ad928SBarry Smith PetscFunctionBegin; 35*4b9ad928SBarry Smith ierr = MatRelax(pc->pmat,x,jac->omega,(MatSORType)flag,0.0,jac->its,jac->lits,y);CHKERRQ(ierr); 36*4b9ad928SBarry Smith PetscFunctionReturn(0); 37*4b9ad928SBarry Smith } 38*4b9ad928SBarry Smith 39*4b9ad928SBarry Smith #undef __FUNCT__ 40*4b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_SOR" 41*4b9ad928SBarry Smith static int PCApplyRichardson_SOR(PC pc,Vec b,Vec y,Vec w,PetscReal rtol,PetscReal atol, PetscReal dtol,int its) 42*4b9ad928SBarry Smith { 43*4b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 44*4b9ad928SBarry Smith int ierr; 45*4b9ad928SBarry Smith 46*4b9ad928SBarry Smith PetscFunctionBegin; 47*4b9ad928SBarry Smith PetscLogInfo(pc,"PCApplyRichardson_SOR: Warning, convergence critera ignored, using %d iterations\n",its); 48*4b9ad928SBarry Smith its = its*jac->its; 49*4b9ad928SBarry Smith ierr = MatRelax(pc->pmat,b,jac->omega,(MatSORType)jac->sym,0.0,its,jac->lits,y);CHKERRQ(ierr); 50*4b9ad928SBarry Smith PetscFunctionReturn(0); 51*4b9ad928SBarry Smith } 52*4b9ad928SBarry Smith 53*4b9ad928SBarry Smith #undef __FUNCT__ 54*4b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_SOR" 55*4b9ad928SBarry Smith static int PCSetFromOptions_SOR(PC pc) 56*4b9ad928SBarry Smith { 57*4b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 58*4b9ad928SBarry Smith int ierr; 59*4b9ad928SBarry Smith PetscTruth flg; 60*4b9ad928SBarry Smith 61*4b9ad928SBarry Smith PetscFunctionBegin; 62*4b9ad928SBarry Smith ierr = PetscOptionsHead("(S)SOR options");CHKERRQ(ierr); 63*4b9ad928SBarry Smith ierr = PetscOptionsReal("-pc_sor_omega","relaxation factor (0 < omega < 2)","PCSORSetOmega",jac->omega,&jac->omega,0);CHKERRQ(ierr); 64*4b9ad928SBarry Smith ierr = PetscOptionsInt("-pc_sor_its","number of inner SOR iterations","PCSORSetIterations",jac->its,&jac->its,0);CHKERRQ(ierr); 65*4b9ad928SBarry Smith ierr = PetscOptionsInt("-pc_sor_lits","number of local inner SOR iterations","PCSORSetIterations",jac->lits,&jac->lits,0);CHKERRQ(ierr); 66*4b9ad928SBarry Smith ierr = PetscOptionsLogicalGroupBegin("-pc_sor_symmetric","SSOR, not SOR","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 67*4b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 68*4b9ad928SBarry Smith ierr = PetscOptionsLogicalGroup("-pc_sor_backward","use backward sweep instead of forward","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 69*4b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_BACKWARD_SWEEP);CHKERRQ(ierr);} 70*4b9ad928SBarry Smith ierr = PetscOptionsLogicalGroup("-pc_sor_local_symmetric","use SSOR seperately on each processor","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 71*4b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_SYMMETRIC_SWEEP);CHKERRQ(ierr);} 72*4b9ad928SBarry Smith ierr = PetscOptionsLogicalGroup("-pc_sor_local_backward","use backward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 73*4b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_BACKWARD_SWEEP);CHKERRQ(ierr);} 74*4b9ad928SBarry Smith ierr = PetscOptionsLogicalGroupEnd("-pc_sor_local_forward","use forward sweep locally","PCSORSetSymmetric",&flg);CHKERRQ(ierr); 75*4b9ad928SBarry Smith if (flg) {ierr = PCSORSetSymmetric(pc,SOR_LOCAL_FORWARD_SWEEP);CHKERRQ(ierr);} 76*4b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 77*4b9ad928SBarry Smith PetscFunctionReturn(0); 78*4b9ad928SBarry Smith } 79*4b9ad928SBarry Smith 80*4b9ad928SBarry Smith #undef __FUNCT__ 81*4b9ad928SBarry Smith #define __FUNCT__ "PCView_SOR" 82*4b9ad928SBarry Smith static int PCView_SOR(PC pc,PetscViewer viewer) 83*4b9ad928SBarry Smith { 84*4b9ad928SBarry Smith PC_SOR *jac = (PC_SOR*)pc->data; 85*4b9ad928SBarry Smith MatSORType sym = jac->sym; 86*4b9ad928SBarry Smith char *sortype; 87*4b9ad928SBarry Smith int ierr; 88*4b9ad928SBarry Smith PetscTruth isascii; 89*4b9ad928SBarry Smith 90*4b9ad928SBarry Smith PetscFunctionBegin; 91*4b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 92*4b9ad928SBarry Smith if (isascii) { 93*4b9ad928SBarry Smith if (sym & SOR_ZERO_INITIAL_GUESS) {ierr = PetscViewerASCIIPrintf(viewer," SOR: zero initial guess\n");CHKERRQ(ierr);} 94*4b9ad928SBarry Smith if (sym == SOR_APPLY_UPPER) sortype = "apply_upper"; 95*4b9ad928SBarry Smith else if (sym == SOR_APPLY_LOWER) sortype = "apply_lower"; 96*4b9ad928SBarry Smith else if (sym & SOR_EISENSTAT) sortype = "Eisenstat"; 97*4b9ad928SBarry Smith else if ((sym & SOR_SYMMETRIC_SWEEP) == SOR_SYMMETRIC_SWEEP) 98*4b9ad928SBarry Smith sortype = "symmetric"; 99*4b9ad928SBarry Smith else if (sym & SOR_BACKWARD_SWEEP) sortype = "backward"; 100*4b9ad928SBarry Smith else if (sym & SOR_FORWARD_SWEEP) sortype = "forward"; 101*4b9ad928SBarry Smith else if ((sym & SOR_LOCAL_SYMMETRIC_SWEEP) == SOR_LOCAL_SYMMETRIC_SWEEP) 102*4b9ad928SBarry Smith sortype = "local_symmetric"; 103*4b9ad928SBarry Smith else if (sym & SOR_LOCAL_FORWARD_SWEEP) sortype = "local_forward"; 104*4b9ad928SBarry Smith else if (sym & SOR_LOCAL_BACKWARD_SWEEP) sortype = "local_backward"; 105*4b9ad928SBarry Smith else sortype = "unknown"; 106*4b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," SOR: type = %s, iterations = %d, omega = %g\n",sortype,jac->its,jac->omega);CHKERRQ(ierr); 107*4b9ad928SBarry Smith } else { 108*4b9ad928SBarry Smith SETERRQ1(1,"Viewer type %s not supported for PCSOR",((PetscObject)viewer)->type_name); 109*4b9ad928SBarry Smith } 110*4b9ad928SBarry Smith PetscFunctionReturn(0); 111*4b9ad928SBarry Smith } 112*4b9ad928SBarry Smith 113*4b9ad928SBarry Smith 114*4b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 115*4b9ad928SBarry Smith EXTERN_C_BEGIN 116*4b9ad928SBarry Smith #undef __FUNCT__ 117*4b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric_SOR" 118*4b9ad928SBarry Smith int PCSORSetSymmetric_SOR(PC pc,MatSORType flag) 119*4b9ad928SBarry Smith { 120*4b9ad928SBarry Smith PC_SOR *jac; 121*4b9ad928SBarry Smith 122*4b9ad928SBarry Smith PetscFunctionBegin; 123*4b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 124*4b9ad928SBarry Smith jac->sym = flag; 125*4b9ad928SBarry Smith PetscFunctionReturn(0); 126*4b9ad928SBarry Smith } 127*4b9ad928SBarry Smith EXTERN_C_END 128*4b9ad928SBarry Smith 129*4b9ad928SBarry Smith EXTERN_C_BEGIN 130*4b9ad928SBarry Smith #undef __FUNCT__ 131*4b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega_SOR" 132*4b9ad928SBarry Smith int PCSORSetOmega_SOR(PC pc,PetscReal omega) 133*4b9ad928SBarry Smith { 134*4b9ad928SBarry Smith PC_SOR *jac; 135*4b9ad928SBarry Smith 136*4b9ad928SBarry Smith PetscFunctionBegin; 137*4b9ad928SBarry Smith if (omega >= 2.0 || omega <= 0.0) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Relaxation out of range"); 138*4b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 139*4b9ad928SBarry Smith jac->omega = omega; 140*4b9ad928SBarry Smith PetscFunctionReturn(0); 141*4b9ad928SBarry Smith } 142*4b9ad928SBarry Smith EXTERN_C_END 143*4b9ad928SBarry Smith 144*4b9ad928SBarry Smith EXTERN_C_BEGIN 145*4b9ad928SBarry Smith #undef __FUNCT__ 146*4b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations_SOR" 147*4b9ad928SBarry Smith int PCSORSetIterations_SOR(PC pc,int its,int lits) 148*4b9ad928SBarry Smith { 149*4b9ad928SBarry Smith PC_SOR *jac; 150*4b9ad928SBarry Smith 151*4b9ad928SBarry Smith PetscFunctionBegin; 152*4b9ad928SBarry Smith jac = (PC_SOR*)pc->data; 153*4b9ad928SBarry Smith jac->its = its; 154*4b9ad928SBarry Smith jac->lits = lits; 155*4b9ad928SBarry Smith PetscFunctionReturn(0); 156*4b9ad928SBarry Smith } 157*4b9ad928SBarry Smith EXTERN_C_END 158*4b9ad928SBarry Smith 159*4b9ad928SBarry Smith /* ------------------------------------------------------------------------------*/ 160*4b9ad928SBarry Smith #undef __FUNCT__ 161*4b9ad928SBarry Smith #define __FUNCT__ "PCSORSetSymmetric" 162*4b9ad928SBarry Smith /*@ 163*4b9ad928SBarry Smith PCSORSetSymmetric - Sets the SOR preconditioner to use symmetric (SSOR), 164*4b9ad928SBarry Smith backward, or forward relaxation. The local variants perform SOR on 165*4b9ad928SBarry Smith each processor. By default forward relaxation is used. 166*4b9ad928SBarry Smith 167*4b9ad928SBarry Smith Collective on PC 168*4b9ad928SBarry Smith 169*4b9ad928SBarry Smith Input Parameters: 170*4b9ad928SBarry Smith + pc - the preconditioner context 171*4b9ad928SBarry Smith - flag - one of the following 172*4b9ad928SBarry Smith .vb 173*4b9ad928SBarry Smith SOR_FORWARD_SWEEP 174*4b9ad928SBarry Smith SOR_BACKWARD_SWEEP 175*4b9ad928SBarry Smith SOR_SYMMETRIC_SWEEP 176*4b9ad928SBarry Smith SOR_LOCAL_FORWARD_SWEEP 177*4b9ad928SBarry Smith SOR_LOCAL_BACKWARD_SWEEP 178*4b9ad928SBarry Smith SOR_LOCAL_SYMMETRIC_SWEEP 179*4b9ad928SBarry Smith .ve 180*4b9ad928SBarry Smith 181*4b9ad928SBarry Smith Options Database Keys: 182*4b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 183*4b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 184*4b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 185*4b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 186*4b9ad928SBarry Smith - -pc_sor_local_backward - Activates local backward version 187*4b9ad928SBarry Smith 188*4b9ad928SBarry Smith Notes: 189*4b9ad928SBarry Smith To use the Eisenstat trick with SSOR, employ the PCEISENSTAT preconditioner, 190*4b9ad928SBarry Smith which can be chosen with the option 191*4b9ad928SBarry Smith . -pc_type eisenstat - Activates Eisenstat trick 192*4b9ad928SBarry Smith 193*4b9ad928SBarry Smith Level: intermediate 194*4b9ad928SBarry Smith 195*4b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, sweep, forward, backward, symmetric 196*4b9ad928SBarry Smith 197*4b9ad928SBarry Smith .seealso: PCEisenstatSetOmega(), PCSORSetIterations(), PCSORSetOmega() 198*4b9ad928SBarry Smith @*/ 199*4b9ad928SBarry Smith int PCSORSetSymmetric(PC pc,MatSORType flag) 200*4b9ad928SBarry Smith { 201*4b9ad928SBarry Smith int ierr,(*f)(PC,MatSORType); 202*4b9ad928SBarry Smith 203*4b9ad928SBarry Smith PetscFunctionBegin; 204*4b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 205*4b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetSymmetric_C",(void (**)(void))&f);CHKERRQ(ierr); 206*4b9ad928SBarry Smith if (f) { 207*4b9ad928SBarry Smith ierr = (*f)(pc,flag);CHKERRQ(ierr); 208*4b9ad928SBarry Smith } 209*4b9ad928SBarry Smith PetscFunctionReturn(0); 210*4b9ad928SBarry Smith } 211*4b9ad928SBarry Smith 212*4b9ad928SBarry Smith #undef __FUNCT__ 213*4b9ad928SBarry Smith #define __FUNCT__ "PCSORSetOmega" 214*4b9ad928SBarry Smith /*@ 215*4b9ad928SBarry Smith PCSORSetOmega - Sets the SOR relaxation coefficient, omega 216*4b9ad928SBarry Smith (where omega = 1.0 by default). 217*4b9ad928SBarry Smith 218*4b9ad928SBarry Smith Collective on PC 219*4b9ad928SBarry Smith 220*4b9ad928SBarry Smith Input Parameters: 221*4b9ad928SBarry Smith + pc - the preconditioner context 222*4b9ad928SBarry Smith - omega - relaxation coefficient (0 < omega < 2). 223*4b9ad928SBarry Smith 224*4b9ad928SBarry Smith Options Database Key: 225*4b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 226*4b9ad928SBarry Smith 227*4b9ad928SBarry Smith Level: intermediate 228*4b9ad928SBarry Smith 229*4b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, relaxation, omega 230*4b9ad928SBarry Smith 231*4b9ad928SBarry Smith .seealso: PCSORSetSymmetric(), PCSORSetIterations(), PCEisenstatSetOmega() 232*4b9ad928SBarry Smith @*/ 233*4b9ad928SBarry Smith int PCSORSetOmega(PC pc,PetscReal omega) 234*4b9ad928SBarry Smith { 235*4b9ad928SBarry Smith int ierr,(*f)(PC,PetscReal); 236*4b9ad928SBarry Smith 237*4b9ad928SBarry Smith PetscFunctionBegin; 238*4b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetOmega_C",(void (**)(void))&f);CHKERRQ(ierr); 239*4b9ad928SBarry Smith if (f) { 240*4b9ad928SBarry Smith ierr = (*f)(pc,omega);CHKERRQ(ierr); 241*4b9ad928SBarry Smith } 242*4b9ad928SBarry Smith PetscFunctionReturn(0); 243*4b9ad928SBarry Smith } 244*4b9ad928SBarry Smith 245*4b9ad928SBarry Smith #undef __FUNCT__ 246*4b9ad928SBarry Smith #define __FUNCT__ "PCSORSetIterations" 247*4b9ad928SBarry Smith /*@ 248*4b9ad928SBarry Smith PCSORSetIterations - Sets the number of inner iterations to 249*4b9ad928SBarry Smith be used by the SOR preconditioner. The default is 1. 250*4b9ad928SBarry Smith 251*4b9ad928SBarry Smith Collective on PC 252*4b9ad928SBarry Smith 253*4b9ad928SBarry Smith Input Parameters: 254*4b9ad928SBarry Smith + pc - the preconditioner context 255*4b9ad928SBarry Smith . lits - number of local iterations, smoothings over just variables on processor 256*4b9ad928SBarry Smith - its - number of parallel iterations to use; each parallel iteration has lits local iterations 257*4b9ad928SBarry Smith 258*4b9ad928SBarry Smith Options Database Key: 259*4b9ad928SBarry Smith + -pc_sor_its <its> - Sets number of iterations 260*4b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 261*4b9ad928SBarry Smith 262*4b9ad928SBarry Smith Level: intermediate 263*4b9ad928SBarry Smith 264*4b9ad928SBarry Smith Notes: When run on one processor the number of smoothings is lits*its 265*4b9ad928SBarry Smith 266*4b9ad928SBarry Smith .keywords: PC, SOR, SSOR, set, iterations 267*4b9ad928SBarry Smith 268*4b9ad928SBarry Smith .seealso: PCSORSetOmega(), PCSORSetSymmetric() 269*4b9ad928SBarry Smith @*/ 270*4b9ad928SBarry Smith int PCSORSetIterations(PC pc,int its,int lits) 271*4b9ad928SBarry Smith { 272*4b9ad928SBarry Smith int ierr,(*f)(PC,int,int); 273*4b9ad928SBarry Smith 274*4b9ad928SBarry Smith PetscFunctionBegin; 275*4b9ad928SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE); 276*4b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCSORSetIterations_C",(void (**)(void))&f);CHKERRQ(ierr); 277*4b9ad928SBarry Smith if (f) { 278*4b9ad928SBarry Smith ierr = (*f)(pc,its,lits);CHKERRQ(ierr); 279*4b9ad928SBarry Smith } 280*4b9ad928SBarry Smith PetscFunctionReturn(0); 281*4b9ad928SBarry Smith } 282*4b9ad928SBarry Smith 283*4b9ad928SBarry Smith /*MC 284*4b9ad928SBarry Smith PCSOR - (S)SOR (successive over relaxation, Gauss-Seidel) preconditioning 285*4b9ad928SBarry Smith 286*4b9ad928SBarry Smith Options Database Keys: 287*4b9ad928SBarry Smith + -pc_sor_symmetric - Activates symmetric version 288*4b9ad928SBarry Smith . -pc_sor_backward - Activates backward version 289*4b9ad928SBarry Smith . -pc_sor_local_forward - Activates local forward version 290*4b9ad928SBarry Smith . -pc_sor_local_symmetric - Activates local symmetric version 291*4b9ad928SBarry Smith . -pc_sor_local_backward - Activates local backward version 292*4b9ad928SBarry Smith . -pc_sor_omega <omega> - Sets omega 293*4b9ad928SBarry Smith . -pc_sor_its <its> - Sets number of iterations 294*4b9ad928SBarry Smith - -pc_sor_lits <lits> - Sets number of local iterations 295*4b9ad928SBarry Smith 296*4b9ad928SBarry Smith Level: beginner 297*4b9ad928SBarry Smith 298*4b9ad928SBarry Smith Concepts: SOR, preconditioners, Gauss-Seidel 299*4b9ad928SBarry Smith 300*4b9ad928SBarry Smith Notes: Only implemented for the AIJ matrix format. 301*4b9ad928SBarry Smith Not a true parallel SOR, in parallel this implementation corresponds to block 302*4b9ad928SBarry Smith Jacobi with SOR on each block. 303*4b9ad928SBarry Smith 304*4b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 305*4b9ad928SBarry Smith PCSORSetIterations(), PCSORSetSymmetric(), PCSORSetOmega(), PCEISENSTAT 306*4b9ad928SBarry Smith M*/ 307*4b9ad928SBarry Smith 308*4b9ad928SBarry Smith EXTERN_C_BEGIN 309*4b9ad928SBarry Smith #undef __FUNCT__ 310*4b9ad928SBarry Smith #define __FUNCT__ "PCCreate_SOR" 311*4b9ad928SBarry Smith int PCCreate_SOR(PC pc) 312*4b9ad928SBarry Smith { 313*4b9ad928SBarry Smith int ierr; 314*4b9ad928SBarry Smith PC_SOR *jac; 315*4b9ad928SBarry Smith 316*4b9ad928SBarry Smith PetscFunctionBegin; 317*4b9ad928SBarry Smith ierr = PetscNew(PC_SOR,&jac);CHKERRQ(ierr); 318*4b9ad928SBarry Smith PetscLogObjectMemory(pc,sizeof(PC_SOR)); 319*4b9ad928SBarry Smith 320*4b9ad928SBarry Smith pc->ops->apply = PCApply_SOR; 321*4b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_SOR; 322*4b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_SOR; 323*4b9ad928SBarry Smith pc->ops->setup = 0; 324*4b9ad928SBarry Smith pc->ops->view = PCView_SOR; 325*4b9ad928SBarry Smith pc->ops->destroy = PCDestroy_SOR; 326*4b9ad928SBarry Smith pc->data = (void*)jac; 327*4b9ad928SBarry Smith jac->sym = SOR_FORWARD_SWEEP; 328*4b9ad928SBarry Smith jac->omega = 1.0; 329*4b9ad928SBarry Smith jac->its = 1; 330*4b9ad928SBarry Smith jac->lits = 1; 331*4b9ad928SBarry Smith 332*4b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetSymmetric_C","PCSORSetSymmetric_SOR", 333*4b9ad928SBarry Smith PCSORSetSymmetric_SOR);CHKERRQ(ierr); 334*4b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetOmega_C","PCSORSetOmega_SOR", 335*4b9ad928SBarry Smith PCSORSetOmega_SOR);CHKERRQ(ierr); 336*4b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCSORSetIterations_C","PCSORSetIterations_SOR", 337*4b9ad928SBarry Smith PCSORSetIterations_SOR);CHKERRQ(ierr); 338*4b9ad928SBarry Smith 339*4b9ad928SBarry Smith PetscFunctionReturn(0); 340*4b9ad928SBarry Smith } 341*4b9ad928SBarry Smith EXTERN_C_END 342*4b9ad928SBarry Smith 343*4b9ad928SBarry Smith 344*4b9ad928SBarry Smith 345*4b9ad928SBarry Smith 346*4b9ad928SBarry Smith 347