1 2 #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscksp.h" I*/ 3 /*I "petscpcmg.h" I*/ 4 5 #undef __FUNCT__ 6 #define __FUNCT__ "PCMGDefaultResidual" 7 /*@C 8 PCMGDefaultResidual - Default routine to calculate the residual. 9 10 Collective on Mat and Vec 11 12 Input Parameters: 13 + mat - the matrix 14 . b - the right-hand-side 15 - x - the approximate solution 16 17 Output Parameter: 18 . r - location to store the residual 19 20 Level: advanced 21 22 .keywords: MG, default, multigrid, residual 23 24 .seealso: PCMGSetResidual() 25 @*/ 26 PetscErrorCode PCMGDefaultResidual(Mat mat,Vec b,Vec x,Vec r) 27 { 28 PetscErrorCode ierr; 29 30 PetscFunctionBegin; 31 ierr = MatMult(mat,x,r);CHKERRQ(ierr); 32 ierr = VecAYPX(r,-1.0,b);CHKERRQ(ierr); 33 PetscFunctionReturn(0); 34 } 35 36 /* ---------------------------------------------------------------------------*/ 37 38 #undef __FUNCT__ 39 #define __FUNCT__ "PCMGGetCoarseSolve" 40 /*@ 41 PCMGGetCoarseSolve - Gets the solver context to be used on the coarse grid. 42 43 Not Collective 44 45 Input Parameter: 46 . pc - the multigrid context 47 48 Output Parameter: 49 . ksp - the coarse grid solver context 50 51 Level: advanced 52 53 .keywords: MG, multigrid, get, coarse grid 54 @*/ 55 PetscErrorCode PCMGGetCoarseSolve(PC pc,KSP *ksp) 56 { 57 PC_MG *mg = (PC_MG*)pc->data; 58 PC_MG_Levels **mglevels = mg->levels; 59 60 PetscFunctionBegin; 61 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 62 *ksp = mglevels[0]->smoothd; 63 PetscFunctionReturn(0); 64 } 65 66 #undef __FUNCT__ 67 #define __FUNCT__ "PCMGSetResidual" 68 /*@C 69 PCMGSetResidual - Sets the function to be used to calculate the residual 70 on the lth level. 71 72 Logically Collective on PC and Mat 73 74 Input Parameters: 75 + pc - the multigrid context 76 . l - the level (0 is coarsest) to supply 77 . residual - function used to form residual (usually PCMGDefaultResidual) 78 - mat - matrix associated with residual 79 80 Level: advanced 81 82 .keywords: MG, set, multigrid, residual, level 83 84 .seealso: PCMGDefaultResidual() 85 @*/ 86 PetscErrorCode PCMGSetResidual(PC pc,PetscInt l,PetscErrorCode (*residual)(Mat,Vec,Vec,Vec),Mat mat) 87 { 88 PC_MG *mg = (PC_MG*)pc->data; 89 PC_MG_Levels **mglevels = mg->levels; 90 PetscErrorCode ierr; 91 92 PetscFunctionBegin; 93 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 94 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 95 mglevels[l]->residual = residual; 96 if (mat) {ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);} 97 ierr = MatDestroy(&mglevels[l]->A);CHKERRQ(ierr); 98 mglevels[l]->A = mat; 99 PetscFunctionReturn(0); 100 } 101 102 #undef __FUNCT__ 103 #define __FUNCT__ "PCMGSetInterpolation" 104 /*@ 105 PCMGSetInterpolation - Sets the function to be used to calculate the 106 interpolation from l-1 to the lth level 107 108 Logically Collective on PC and Mat 109 110 Input Parameters: 111 + pc - the multigrid context 112 . mat - the interpolation operator 113 - l - the level (0 is coarsest) to supply [do not supply 0] 114 115 Level: advanced 116 117 Notes: 118 Usually this is the same matrix used also to set the restriction 119 for the same level. 120 121 One can pass in the interpolation matrix or its transpose; PETSc figures 122 out from the matrix size which one it is. 123 124 .keywords: multigrid, set, interpolate, level 125 126 .seealso: PCMGSetRestriction() 127 @*/ 128 PetscErrorCode PCMGSetInterpolation(PC pc,PetscInt l,Mat mat) 129 { 130 PC_MG *mg = (PC_MG*)pc->data; 131 PC_MG_Levels **mglevels = mg->levels; 132 PetscErrorCode ierr; 133 134 PetscFunctionBegin; 135 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 136 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 137 if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set interpolation routine for coarsest level"); 138 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 139 ierr = MatDestroy(&mglevels[l]->interpolate);CHKERRQ(ierr); 140 mglevels[l]->interpolate = mat; 141 PetscFunctionReturn(0); 142 } 143 144 #undef __FUNCT__ 145 #define __FUNCT__ "PCMGSetRestriction" 146 /*@ 147 PCMGSetRestriction - Sets the function to be used to restrict vector 148 from level l to l-1. 149 150 Logically Collective on PC and Mat 151 152 Input Parameters: 153 + pc - the multigrid context 154 . mat - the restriction matrix 155 - l - the level (0 is coarsest) to supply [Do not supply 0] 156 157 Level: advanced 158 159 Notes: 160 Usually this is the same matrix used also to set the interpolation 161 for the same level. 162 163 One can pass in the interpolation matrix or its transpose; PETSc figures 164 out from the matrix size which one it is. 165 166 If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 167 is used. 168 169 .keywords: MG, set, multigrid, restriction, level 170 171 .seealso: PCMGSetInterpolation() 172 @*/ 173 PetscErrorCode PCMGSetRestriction(PC pc,PetscInt l,Mat mat) 174 { 175 PetscErrorCode ierr; 176 PC_MG *mg = (PC_MG*)pc->data; 177 PC_MG_Levels **mglevels = mg->levels; 178 179 PetscFunctionBegin; 180 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 181 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 182 if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 183 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 184 ierr = MatDestroy(&mglevels[l]->restrct);CHKERRQ(ierr); 185 mglevels[l]->restrct = mat; 186 PetscFunctionReturn(0); 187 } 188 189 #undef __FUNCT__ 190 #define __FUNCT__ "PCMGSetRScale" 191 /*@ 192 PCMGSetRScale - Sets the pointwise scaling for the restriction operator from level l to l-1. 193 194 Logically Collective on PC and Mat 195 196 Input Parameters: 197 + pc - the multigrid context 198 . rscale - the scaling 199 - l - the level (0 is coarsest) to supply [Do not supply 0] 200 201 Level: advanced 202 203 Notes: 204 When evaluating a function on a coarse level one does not want to do F( R * x) one does F( rscale * R * x) where rscale is 1 over the row sums of R. 205 206 .keywords: MG, set, multigrid, restriction, level 207 208 .seealso: PCMGSetInterpolation(), PCMGSetRestriction() 209 @*/ 210 PetscErrorCode PCMGSetRScale(PC pc,PetscInt l,Vec rscale) 211 { 212 PetscErrorCode ierr; 213 PC_MG *mg = (PC_MG*)pc->data; 214 PC_MG_Levels **mglevels = mg->levels; 215 216 PetscFunctionBegin; 217 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 218 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 219 if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Do not set restriction routine for coarsest level"); 220 ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 221 ierr = VecDestroy(&mglevels[l]->rscale);CHKERRQ(ierr); 222 mglevels[l]->rscale = rscale; 223 PetscFunctionReturn(0); 224 } 225 226 #undef __FUNCT__ 227 #define __FUNCT__ "PCMGGetSmoother" 228 /*@ 229 PCMGGetSmoother - Gets the KSP context to be used as smoother for 230 both pre- and post-smoothing. Call both PCMGGetSmootherUp() and 231 PCMGGetSmootherDown() to use different functions for pre- and 232 post-smoothing. 233 234 Not Collective, KSP returned is parallel if PC is 235 236 Input Parameters: 237 + pc - the multigrid context 238 - l - the level (0 is coarsest) to supply 239 240 Ouput Parameters: 241 . ksp - the smoother 242 243 Level: advanced 244 245 .keywords: MG, get, multigrid, level, smoother, pre-smoother, post-smoother 246 247 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 248 @*/ 249 PetscErrorCode PCMGGetSmoother(PC pc,PetscInt l,KSP *ksp) 250 { 251 PC_MG *mg = (PC_MG*)pc->data; 252 PC_MG_Levels **mglevels = mg->levels; 253 254 PetscFunctionBegin; 255 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 256 *ksp = mglevels[l]->smoothd; 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "PCMGGetSmootherUp" 262 /*@ 263 PCMGGetSmootherUp - Gets the KSP context to be used as smoother after 264 coarse grid correction (post-smoother). 265 266 Not Collective, KSP returned is parallel if PC is 267 268 Input Parameters: 269 + pc - the multigrid context 270 - l - the level (0 is coarsest) to supply 271 272 Ouput Parameters: 273 . ksp - the smoother 274 275 Level: advanced 276 277 .keywords: MG, multigrid, get, smoother, up, post-smoother, level 278 279 .seealso: PCMGGetSmootherUp(), PCMGGetSmootherDown() 280 @*/ 281 PetscErrorCode PCMGGetSmootherUp(PC pc,PetscInt l,KSP *ksp) 282 { 283 PC_MG *mg = (PC_MG*)pc->data; 284 PC_MG_Levels **mglevels = mg->levels; 285 PetscErrorCode ierr; 286 const char *prefix; 287 MPI_Comm comm; 288 289 PetscFunctionBegin; 290 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 291 /* 292 This is called only if user wants a different pre-smoother from post. 293 Thus we check if a different one has already been allocated, 294 if not we allocate it. 295 */ 296 if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"There is no such thing as a up smoother on the coarse grid"); 297 if (mglevels[l]->smoothu == mglevels[l]->smoothd) { 298 ierr = PetscObjectGetComm((PetscObject)mglevels[l]->smoothd,&comm);CHKERRQ(ierr); 299 ierr = KSPGetOptionsPrefix(mglevels[l]->smoothd,&prefix);CHKERRQ(ierr); 300 ierr = KSPCreate(comm,&mglevels[l]->smoothu);CHKERRQ(ierr); 301 ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[l]->smoothu,(PetscObject)pc,mglevels[0]->levels-l);CHKERRQ(ierr); 302 ierr = KSPSetTolerances(mglevels[l]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 303 ierr = KSPSetOptionsPrefix(mglevels[l]->smoothu,prefix);CHKERRQ(ierr); 304 ierr = PetscLogObjectParent(pc,mglevels[l]->smoothu);CHKERRQ(ierr); 305 } 306 if (ksp) *ksp = mglevels[l]->smoothu; 307 PetscFunctionReturn(0); 308 } 309 310 #undef __FUNCT__ 311 #define __FUNCT__ "PCMGGetSmootherDown" 312 /*@ 313 PCMGGetSmootherDown - Gets the KSP context to be used as smoother before 314 coarse grid correction (pre-smoother). 315 316 Not Collective, KSP returned is parallel if PC is 317 318 Input Parameters: 319 + pc - the multigrid context 320 - l - the level (0 is coarsest) to supply 321 322 Ouput Parameters: 323 . ksp - the smoother 324 325 Level: advanced 326 327 .keywords: MG, multigrid, get, smoother, down, pre-smoother, level 328 329 .seealso: PCMGGetSmootherUp(), PCMGGetSmoother() 330 @*/ 331 PetscErrorCode PCMGGetSmootherDown(PC pc,PetscInt l,KSP *ksp) 332 { 333 PetscErrorCode ierr; 334 PC_MG *mg = (PC_MG*)pc->data; 335 PC_MG_Levels **mglevels = mg->levels; 336 337 PetscFunctionBegin; 338 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 339 /* make sure smoother up and down are different */ 340 if (l) { 341 ierr = PCMGGetSmootherUp(pc,l,PETSC_NULL);CHKERRQ(ierr); 342 } 343 *ksp = mglevels[l]->smoothd; 344 PetscFunctionReturn(0); 345 } 346 347 #undef __FUNCT__ 348 #define __FUNCT__ "PCMGSetCyclesOnLevel" 349 /*@ 350 PCMGSetCyclesOnLevel - Sets the number of cycles to run on this level. 351 352 Logically Collective on PC 353 354 Input Parameters: 355 + pc - the multigrid context 356 . l - the level (0 is coarsest) this is to be used for 357 - n - the number of cycles 358 359 Level: advanced 360 361 .keywords: MG, multigrid, set, cycles, V-cycle, W-cycle, level 362 363 .seealso: PCMGSetCycles() 364 @*/ 365 PetscErrorCode PCMGSetCyclesOnLevel(PC pc,PetscInt l,PetscInt c) 366 { 367 PC_MG *mg = (PC_MG*)pc->data; 368 PC_MG_Levels **mglevels = mg->levels; 369 370 PetscFunctionBegin; 371 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 372 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 373 PetscValidLogicalCollectiveInt(pc,l,2); 374 PetscValidLogicalCollectiveInt(pc,c,3); 375 mglevels[l]->cycles = c; 376 PetscFunctionReturn(0); 377 } 378 379 #undef __FUNCT__ 380 #define __FUNCT__ "PCMGSetRhs" 381 /*@ 382 PCMGSetRhs - Sets the vector space to be used to store the right-hand side 383 on a particular level. 384 385 Logically Collective on PC and Vec 386 387 Input Parameters: 388 + pc - the multigrid context 389 . l - the level (0 is coarsest) this is to be used for 390 - c - the space 391 392 Level: advanced 393 394 Notes: If this is not provided PETSc will automatically generate one. 395 396 You do not need to keep a reference to this vector if you do 397 not need it PCDestroy() will properly free it. 398 399 .keywords: MG, multigrid, set, right-hand-side, rhs, level 400 401 .seealso: PCMGSetX(), PCMGSetR() 402 @*/ 403 PetscErrorCode PCMGSetRhs(PC pc,PetscInt l,Vec c) 404 { 405 PetscErrorCode ierr; 406 PC_MG *mg = (PC_MG*)pc->data; 407 PC_MG_Levels **mglevels = mg->levels; 408 409 PetscFunctionBegin; 410 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 411 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 412 if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set rhs for finest level"); 413 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 414 ierr = VecDestroy(&mglevels[l]->b);CHKERRQ(ierr); 415 mglevels[l]->b = c; 416 PetscFunctionReturn(0); 417 } 418 419 #undef __FUNCT__ 420 #define __FUNCT__ "PCMGSetX" 421 /*@ 422 PCMGSetX - Sets the vector space to be used to store the solution on a 423 particular level. 424 425 Logically Collective on PC and Vec 426 427 Input Parameters: 428 + pc - the multigrid context 429 . l - the level (0 is coarsest) this is to be used for 430 - c - the space 431 432 Level: advanced 433 434 Notes: If this is not provided PETSc will automatically generate one. 435 436 You do not need to keep a reference to this vector if you do 437 not need it PCDestroy() will properly free it. 438 439 .keywords: MG, multigrid, set, solution, level 440 441 .seealso: PCMGSetRhs(), PCMGSetR() 442 @*/ 443 PetscErrorCode PCMGSetX(PC pc,PetscInt l,Vec c) 444 { 445 PetscErrorCode ierr; 446 PC_MG *mg = (PC_MG*)pc->data; 447 PC_MG_Levels **mglevels = mg->levels; 448 449 PetscFunctionBegin; 450 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 451 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 452 if (l == mglevels[0]->levels-1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_INCOMP,"Do not set x for finest level"); 453 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 454 ierr = VecDestroy(&mglevels[l]->x);CHKERRQ(ierr); 455 mglevels[l]->x = c; 456 PetscFunctionReturn(0); 457 } 458 459 #undef __FUNCT__ 460 #define __FUNCT__ "PCMGSetR" 461 /*@ 462 PCMGSetR - Sets the vector space to be used to store the residual on a 463 particular level. 464 465 Logically Collective on PC and Vec 466 467 Input Parameters: 468 + pc - the multigrid context 469 . l - the level (0 is coarsest) this is to be used for 470 - c - the space 471 472 Level: advanced 473 474 Notes: If this is not provided PETSc will automatically generate one. 475 476 You do not need to keep a reference to this vector if you do 477 not need it PCDestroy() will properly free it. 478 479 .keywords: MG, multigrid, set, residual, level 480 @*/ 481 PetscErrorCode PCMGSetR(PC pc,PetscInt l,Vec c) 482 { 483 PetscErrorCode ierr; 484 PC_MG *mg = (PC_MG*)pc->data; 485 PC_MG_Levels **mglevels = mg->levels; 486 487 PetscFunctionBegin; 488 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 489 if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 490 if (!l) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Need not set residual vector for coarse grid"); 491 ierr = PetscObjectReference((PetscObject)c);CHKERRQ(ierr); 492 ierr = VecDestroy(&mglevels[l]->r);CHKERRQ(ierr); 493 mglevels[l]->r = c; 494 PetscFunctionReturn(0); 495 } 496