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