1 /* 2 Defines the multigrid preconditioner interface. 3 */ 4 #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 5 6 7 /* 8 MGMCycle_Private - Given an MG structure created with MGCreate() runs 9 one multiplicative cycle down through the levels and 10 back up. 11 12 Input Parameter: 13 . mg - structure created with MGCreate(). 14 */ 15 #undef __FUNCT__ 16 #define __FUNCT__ "MGMCycle_Private" 17 PetscErrorCode MGMCycle_Private(MG *mglevels,PetscTruth *converged) 18 { 19 MG mg = *mglevels,mgc; 20 int cycles = mg->cycles,ierr; 21 PetscScalar zero = 0.0; 22 23 PetscFunctionBegin; 24 if (converged) *converged = PETSC_FALSE; 25 26 if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);} 27 ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr); 28 if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);} 29 if (mg->level) { /* not the coarsest grid */ 30 ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr); 31 32 /* if on finest level and have convergence criteria set */ 33 if (mg->level == mg->levels-1 && mg->ttol) { 34 PetscReal rnorm; 35 ierr = VecNorm(mg->r,NORM_2,&rnorm);CHKERRQ(ierr); 36 if (rnorm <= mg->ttol) { 37 *converged = PETSC_TRUE; 38 if (rnorm < mg->atol) { 39 PetscLogInfo(0,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",rnorm,mg->atol); 40 } else { 41 PetscLogInfo(0,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",rnorm,mg->ttol); 42 } 43 PetscFunctionReturn(0); 44 } 45 } 46 47 mgc = *(mglevels - 1); 48 ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr); 49 ierr = VecSet(&zero,mgc->x);CHKERRQ(ierr); 50 while (cycles--) { 51 ierr = MGMCycle_Private(mglevels-1,converged);CHKERRQ(ierr); 52 } 53 ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr); 54 if (mg->eventsolve) {ierr = PetscLogEventBegin(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);} 55 ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr); 56 if (mg->eventsolve) {ierr = PetscLogEventEnd(mg->eventsolve,0,0,0,0);CHKERRQ(ierr);} 57 } 58 PetscFunctionReturn(0); 59 } 60 61 /* 62 MGCreate_Private - Creates a MG structure for use with the 63 multigrid code. Level 0 is the coarsest. (But the 64 finest level is stored first in the array). 65 66 */ 67 #undef __FUNCT__ 68 #define __FUNCT__ "MGCreate_Private" 69 static int MGCreate_Private(MPI_Comm comm,int levels,PC pc,MPI_Comm *comms,MG **result) 70 { 71 MG *mg; 72 int i,ierr,size; 73 char *prefix; 74 PC ipc; 75 76 PetscFunctionBegin; 77 ierr = PetscMalloc(levels*sizeof(MG),&mg);CHKERRQ(ierr); 78 PetscLogObjectMemory(pc,levels*(sizeof(MG)+sizeof(struct _MG))); 79 80 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 81 82 for (i=0; i<levels; i++) { 83 ierr = PetscNew(struct _MG,&mg[i]);CHKERRQ(ierr); 84 ierr = PetscMemzero(mg[i],sizeof(struct _MG));CHKERRQ(ierr); 85 mg[i]->level = i; 86 mg[i]->levels = levels; 87 mg[i]->cycles = 1; 88 89 if (comms) comm = comms[i]; 90 ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr); 91 ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1);CHKERRQ(ierr); 92 ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr); 93 94 /* do special stuff for coarse grid */ 95 if (!i && levels > 1) { 96 ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 97 98 /* coarse solve is (redundant) LU by default */ 99 ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 100 ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr); 101 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 102 if (size > 1) { 103 ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 104 ierr = PCRedundantGetPC(ipc,&ipc);CHKERRQ(ierr); 105 } 106 ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 107 108 } else { 109 ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,"mg_levels_");CHKERRQ(ierr); 110 } 111 PetscLogObjectParent(pc,mg[i]->smoothd); 112 mg[i]->smoothu = mg[i]->smoothd; 113 mg[i]->default_smoothu = 10000; 114 mg[i]->default_smoothd = 10000; 115 mg[i]->rtol = 0.0; 116 mg[i]->atol = 0.0; 117 mg[i]->dtol = 0.0; 118 mg[i]->ttol = 0.0; 119 mg[i]->eventsetup = 0; 120 mg[i]->eventsolve = 0; 121 } 122 *result = mg; 123 PetscFunctionReturn(0); 124 } 125 126 #undef __FUNCT__ 127 #define __FUNCT__ "PCDestroy_MG" 128 static int PCDestroy_MG(PC pc) 129 { 130 MG *mg = (MG*)pc->data; 131 int i,n = mg[0]->levels,ierr; 132 133 PetscFunctionBegin; 134 for (i=0; i<n; i++) { 135 if (mg[i]->smoothd != mg[i]->smoothu) { 136 ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr); 137 } 138 ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr); 139 ierr = PetscFree(mg[i]);CHKERRQ(ierr); 140 } 141 ierr = PetscFree(mg);CHKERRQ(ierr); 142 PetscFunctionReturn(0); 143 } 144 145 146 147 EXTERN PetscErrorCode MGACycle_Private(MG*); 148 EXTERN PetscErrorCode MGFCycle_Private(MG*); 149 EXTERN PetscErrorCode MGKCycle_Private(MG*); 150 151 /* 152 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 153 or full cycle of multigrid. 154 155 Note: 156 A simple wrapper which calls MGMCycle(),MGACycle(), or MGFCycle(). 157 */ 158 #undef __FUNCT__ 159 #define __FUNCT__ "PCApply_MG" 160 static int PCApply_MG(PC pc,Vec b,Vec x) 161 { 162 MG *mg = (MG*)pc->data; 163 PetscScalar zero = 0.0; 164 int levels = mg[0]->levels,ierr; 165 166 PetscFunctionBegin; 167 mg[levels-1]->b = b; 168 mg[levels-1]->x = x; 169 if (mg[0]->am == MGMULTIPLICATIVE) { 170 ierr = VecSet(&zero,x);CHKERRQ(ierr); 171 ierr = MGMCycle_Private(mg+levels-1,PETSC_NULL);CHKERRQ(ierr); 172 } 173 else if (mg[0]->am == MGADDITIVE) { 174 ierr = MGACycle_Private(mg);CHKERRQ(ierr); 175 } 176 else if (mg[0]->am == MGKASKADE) { 177 ierr = MGKCycle_Private(mg);CHKERRQ(ierr); 178 } 179 else { 180 ierr = MGFCycle_Private(mg);CHKERRQ(ierr); 181 } 182 PetscFunctionReturn(0); 183 } 184 185 #undef __FUNCT__ 186 #define __FUNCT__ "PCApplyRichardson_MG" 187 static int PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal atol, PetscReal dtol,int its) 188 { 189 MG *mg = (MG*)pc->data; 190 PetscErrorCode ierr; 191 int levels = mg[0]->levels; 192 PetscTruth converged = PETSC_FALSE; 193 194 PetscFunctionBegin; 195 mg[levels-1]->b = b; 196 mg[levels-1]->x = x; 197 198 mg[levels-1]->rtol = rtol; 199 mg[levels-1]->atol = atol; 200 mg[levels-1]->dtol = dtol; 201 if (rtol) { 202 /* compute initial residual norm for relative convergence test */ 203 PetscReal rnorm; 204 ierr = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr); 205 ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 206 mg[levels-1]->ttol = PetscMax(rtol*rnorm,atol); 207 } else if (atol) { 208 mg[levels-1]->ttol = atol; 209 } else { 210 mg[levels-1]->ttol = 0.0; 211 } 212 213 while (its-- && !converged) { 214 ierr = MGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr); 215 } 216 PetscFunctionReturn(0); 217 } 218 219 #undef __FUNCT__ 220 #define __FUNCT__ "PCSetFromOptions_MG" 221 static int PCSetFromOptions_MG(PC pc) 222 { 223 PetscErrorCode ierr; 224 int indx,m,levels = 1; 225 PetscTruth flg; 226 const char *type[] = {"additive","multiplicative","full","cascade","kascade"}; 227 228 PetscFunctionBegin; 229 230 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 231 if (!pc->data) { 232 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","MGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 233 ierr = MGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 234 } 235 ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 236 if (flg) { 237 ierr = MGSetCycles(pc,m);CHKERRQ(ierr); 238 } 239 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 240 if (flg) { 241 ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 242 } 243 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 244 if (flg) { 245 ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 246 } 247 ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr); 248 if (flg) { 249 MGType mg = (MGType) 0; 250 switch (indx) { 251 case 0: 252 mg = MGADDITIVE; 253 break; 254 case 1: 255 mg = MGMULTIPLICATIVE; 256 break; 257 case 2: 258 mg = MGFULL; 259 break; 260 case 3: 261 mg = MGKASKADE; 262 break; 263 case 4: 264 mg = MGKASKADE; 265 break; 266 } 267 ierr = MGSetType(pc,mg);CHKERRQ(ierr); 268 } 269 ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr); 270 if (flg) { 271 MG *mg = (MG*)pc->data; 272 int i; 273 char eventname[128]; 274 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 275 levels = mg[0]->levels; 276 for (i=0; i<levels; i++) { 277 sprintf(eventname,"MSetup Level %d",i); 278 ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr); 279 sprintf(eventname,"MGSolve Level %d",i); 280 ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr); 281 } 282 } 283 ierr = PetscOptionsTail();CHKERRQ(ierr); 284 PetscFunctionReturn(0); 285 } 286 287 #undef __FUNCT__ 288 #define __FUNCT__ "PCView_MG" 289 static int PCView_MG(PC pc,PetscViewer viewer) 290 { 291 MG *mg = (MG*)pc->data; 292 PetscErrorCode ierr; 293 int levels = mg[0]->levels,i; 294 const char *cstring; 295 PetscTruth iascii; 296 297 PetscFunctionBegin; 298 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 299 if (iascii) { 300 if (mg[0]->am == MGMULTIPLICATIVE) cstring = "multiplicative"; 301 else if (mg[0]->am == MGADDITIVE) cstring = "additive"; 302 else if (mg[0]->am == MGFULL) cstring = "full"; 303 else if (mg[0]->am == MGKASKADE) cstring = "Kaskade"; 304 else cstring = "unknown"; 305 ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%d cycles=%d, pre-smooths=%d, post-smooths=%d\n", 306 cstring,levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr); 307 for (i=0; i<levels; i++) { 308 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %d -------------------------------\n",i);CHKERRQ(ierr); 309 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 310 ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr); 311 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 312 if (mg[i]->smoothd == mg[i]->smoothu) { 313 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 314 } else { 315 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %d -------------------------------\n",i);CHKERRQ(ierr); 316 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 317 ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr); 318 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 319 } 320 } 321 } else { 322 SETERRQ1(1,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 323 } 324 PetscFunctionReturn(0); 325 } 326 327 /* 328 Calls setup for the KSP on each level 329 */ 330 #undef __FUNCT__ 331 #define __FUNCT__ "PCSetUp_MG" 332 static int PCSetUp_MG(PC pc) 333 { 334 MG *mg = (MG*)pc->data; 335 PetscErrorCode ierr; 336 int i,n = mg[0]->levels; 337 PC cpc; 338 PetscTruth preonly,lu,redundant,monitor = PETSC_FALSE,dump; 339 PetscViewer ascii; 340 MPI_Comm comm; 341 342 PetscFunctionBegin; 343 if (!pc->setupcalled) { 344 ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr); 345 346 for (i=1; i<n; i++) { 347 if (mg[i]->smoothd) { 348 if (monitor) { 349 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr); 350 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 351 ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr); 352 ierr = KSPSetMonitor(mg[i]->smoothd,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 353 } 354 ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr); 355 } 356 } 357 for (i=0; i<n; i++) { 358 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 359 if (monitor) { 360 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr); 361 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 362 ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr); 363 ierr = KSPSetMonitor(mg[i]->smoothu,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 364 } 365 ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr); 366 } 367 } 368 } 369 370 for (i=1; i<n; i++) { 371 if (mg[i]->smoothd) { 372 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 373 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 374 ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 375 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 376 } 377 } 378 for (i=0; i<n; i++) { 379 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 380 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 381 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 382 ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 383 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 384 } 385 } 386 387 /* 388 If coarse solver is not direct method then DO NOT USE preonly 389 */ 390 ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 391 if (preonly) { 392 ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 393 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 394 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 395 if (!lu && !redundant) { 396 ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 397 } 398 } 399 400 if (!pc->setupcalled) { 401 if (monitor) { 402 ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 403 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 404 ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr); 405 ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 406 } 407 ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr); 408 } 409 410 if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 411 ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr); 412 if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 413 414 /* 415 Dump the interpolation/restriction matrices to matlab plus the 416 Jacobian/stiffness on each level. This allows Matlab users to 417 easily check if the Galerkin condition A_c = R A_f R^T is satisfied */ 418 ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr); 419 if (dump) { 420 for (i=1; i<n; i++) { 421 ierr = MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr); 422 } 423 for (i=0; i<n; i++) { 424 ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 425 ierr = MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr); 426 } 427 } 428 429 PetscFunctionReturn(0); 430 } 431 432 /* -------------------------------------------------------------------------------------*/ 433 434 #undef __FUNCT__ 435 #define __FUNCT__ "MGSetLevels" 436 /*@C 437 MGSetLevels - Sets the number of levels to use with MG. 438 Must be called before any other MG routine. 439 440 Collective on PC 441 442 Input Parameters: 443 + pc - the preconditioner context 444 . levels - the number of levels 445 - comms - optional communicators for each level; this is to allow solving the coarser problems 446 on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 447 448 Level: intermediate 449 450 Notes: 451 If the number of levels is one then the multigrid uses the -mg_levels prefix 452 for setting the level options rather than the -mg_coarse prefix. 453 454 .keywords: MG, set, levels, multigrid 455 456 .seealso: MGSetType(), MGGetLevels() 457 @*/ 458 PetscErrorCode MGSetLevels(PC pc,int levels,MPI_Comm *comms) 459 { 460 PetscErrorCode ierr; 461 MG *mg; 462 463 PetscFunctionBegin; 464 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 465 466 if (pc->data) { 467 SETERRQ(1,"Number levels already set for MG\n\ 468 make sure that you call MGSetLevels() before KSPSetFromOptions()"); 469 } 470 ierr = MGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr); 471 mg[0]->am = MGMULTIPLICATIVE; 472 pc->data = (void*)mg; 473 pc->ops->applyrichardson = PCApplyRichardson_MG; 474 PetscFunctionReturn(0); 475 } 476 477 #undef __FUNCT__ 478 #define __FUNCT__ "MGGetLevels" 479 /*@ 480 MGGetLevels - Gets the number of levels to use with MG. 481 482 Not Collective 483 484 Input Parameter: 485 . pc - the preconditioner context 486 487 Output parameter: 488 . levels - the number of levels 489 490 Level: advanced 491 492 .keywords: MG, get, levels, multigrid 493 494 .seealso: MGSetLevels() 495 @*/ 496 PetscErrorCode MGGetLevels(PC pc,int *levels) 497 { 498 MG *mg; 499 500 PetscFunctionBegin; 501 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 502 PetscValidIntPointer(levels,2); 503 504 mg = (MG*)pc->data; 505 *levels = mg[0]->levels; 506 PetscFunctionReturn(0); 507 } 508 509 #undef __FUNCT__ 510 #define __FUNCT__ "MGSetType" 511 /*@ 512 MGSetType - Determines the form of multigrid to use: 513 multiplicative, additive, full, or the Kaskade algorithm. 514 515 Collective on PC 516 517 Input Parameters: 518 + pc - the preconditioner context 519 - form - multigrid form, one of MGMULTIPLICATIVE, MGADDITIVE, 520 MGFULL, MGKASKADE 521 522 Options Database Key: 523 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 524 additive, full, kaskade 525 526 Level: advanced 527 528 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 529 530 .seealso: MGSetLevels() 531 @*/ 532 PetscErrorCode MGSetType(PC pc,MGType form) 533 { 534 MG *mg; 535 536 PetscFunctionBegin; 537 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 538 mg = (MG*)pc->data; 539 540 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 541 mg[0]->am = form; 542 if (form == MGMULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 543 else pc->ops->applyrichardson = 0; 544 PetscFunctionReturn(0); 545 } 546 547 #undef __FUNCT__ 548 #define __FUNCT__ "MGSetCycles" 549 /*@ 550 MGSetCycles - Sets the type cycles to use. Use MGSetCyclesOnLevel() for more 551 complicated cycling. 552 553 Collective on PC 554 555 Input Parameters: 556 + mg - the multigrid context 557 - n - the number of cycles 558 559 Options Database Key: 560 $ -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle. 561 562 Level: advanced 563 564 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 565 566 .seealso: MGSetCyclesOnLevel() 567 @*/ 568 PetscErrorCode MGSetCycles(PC pc,int n) 569 { 570 MG *mg; 571 int i,levels; 572 573 PetscFunctionBegin; 574 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 575 mg = (MG*)pc->data; 576 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 577 levels = mg[0]->levels; 578 579 for (i=0; i<levels; i++) { 580 mg[i]->cycles = n; 581 } 582 PetscFunctionReturn(0); 583 } 584 585 #undef __FUNCT__ 586 #define __FUNCT__ "MGCheck" 587 /*@ 588 MGCheck - Checks that all components of the MG structure have 589 been set. 590 591 Collective on PC 592 593 Input Parameters: 594 . mg - the MG structure 595 596 Level: advanced 597 598 .keywords: MG, check, set, multigrid 599 @*/ 600 PetscErrorCode MGCheck(PC pc) 601 { 602 MG *mg; 603 int i,n,count = 0; 604 605 PetscFunctionBegin; 606 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 607 mg = (MG*)pc->data; 608 609 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 610 611 n = mg[0]->levels; 612 613 for (i=1; i<n; i++) { 614 if (!mg[i]->restrct) { 615 (*PetscErrorPrintf)("No restrict set level %d \n",n-i); count++; 616 } 617 if (!mg[i]->interpolate) { 618 (*PetscErrorPrintf)("No interpolate set level %d \n",n-i); count++; 619 } 620 if (!mg[i]->residual) { 621 (*PetscErrorPrintf)("No residual set level %d \n",n-i); count++; 622 } 623 if (!mg[i]->smoothu) { 624 (*PetscErrorPrintf)("No smoothup set level %d \n",n-i); count++; 625 } 626 if (!mg[i]->smoothd) { 627 (*PetscErrorPrintf)("No smoothdown set level %d \n",n-i); count++; 628 } 629 if (!mg[i]->r) { 630 (*PetscErrorPrintf)("No r set level %d \n",n-i); count++; 631 } 632 if (!mg[i-1]->x) { 633 (*PetscErrorPrintf)("No x set level %d \n",n-i); count++; 634 } 635 if (!mg[i-1]->b) { 636 (*PetscErrorPrintf)("No b set level %d \n",n-i); count++; 637 } 638 } 639 PetscFunctionReturn(count); 640 } 641 642 643 #undef __FUNCT__ 644 #define __FUNCT__ "MGSetNumberSmoothDown" 645 /*@ 646 MGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 647 use on all levels. Use MGGetSmootherDown() to set different 648 pre-smoothing steps on different levels. 649 650 Collective on PC 651 652 Input Parameters: 653 + mg - the multigrid context 654 - n - the number of smoothing steps 655 656 Options Database Key: 657 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 658 659 Level: advanced 660 661 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 662 663 .seealso: MGSetNumberSmoothUp() 664 @*/ 665 PetscErrorCode MGSetNumberSmoothDown(PC pc,int n) 666 { 667 MG *mg; 668 int i,levels,ierr; 669 670 PetscFunctionBegin; 671 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 672 mg = (MG*)pc->data; 673 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 674 levels = mg[0]->levels; 675 676 for (i=0; i<levels; i++) { 677 /* make sure smoother up and down are different */ 678 ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 679 ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 680 mg[i]->default_smoothd = n; 681 } 682 PetscFunctionReturn(0); 683 } 684 685 #undef __FUNCT__ 686 #define __FUNCT__ "MGSetNumberSmoothUp" 687 /*@ 688 MGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 689 on all levels. Use MGGetSmootherUp() to set different numbers of 690 post-smoothing steps on different levels. 691 692 Collective on PC 693 694 Input Parameters: 695 + mg - the multigrid context 696 - n - the number of smoothing steps 697 698 Options Database Key: 699 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 700 701 Level: advanced 702 703 Note: this does not set a value on the coarsest grid, since we assume that 704 there is no seperate smooth up on the coarsest grid. If you want to have a 705 seperate smooth up on the coarsest grid then call MGGetSmoothUp(pc,0,&ksp); 706 707 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 708 709 .seealso: MGSetNumberSmoothDown() 710 @*/ 711 PetscErrorCode MGSetNumberSmoothUp(PC pc,int n) 712 { 713 MG *mg; 714 int i,levels,ierr; 715 716 PetscFunctionBegin; 717 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 718 mg = (MG*)pc->data; 719 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 720 levels = mg[0]->levels; 721 722 for (i=1; i<levels; i++) { 723 /* make sure smoother up and down are different */ 724 ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 725 ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 726 mg[i]->default_smoothu = n; 727 } 728 PetscFunctionReturn(0); 729 } 730 731 /* ----------------------------------------------------------------------------------------*/ 732 733 /*MC 734 PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional 735 information about the coarser grid matrices and restriction/interpolation operators. 736 737 Options Database Keys: 738 + -pc_mg_levels <nlevels> - number of levels including finest 739 . -pc_mg_cycles 1 or 2 - for V or W-cycle 740 . -pc_mg_smoothup <n> - number of smoothing steps before interpolation 741 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 742 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 743 . -pc_mg_log - log information about time spent on each level of the solver 744 . -pc_mg_monitor - print information on the multigrid convergence 745 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 746 to the Socket viewer for reading from Matlab. 747 748 Notes: 749 750 Level: intermediate 751 752 Concepts: multigrid 753 754 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 755 MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(), 756 MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(), 757 MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(), 758 MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR() 759 M*/ 760 761 EXTERN_C_BEGIN 762 #undef __FUNCT__ 763 #define __FUNCT__ "PCCreate_MG" 764 PetscErrorCode PCCreate_MG(PC pc) 765 { 766 PetscFunctionBegin; 767 pc->ops->apply = PCApply_MG; 768 pc->ops->setup = PCSetUp_MG; 769 pc->ops->destroy = PCDestroy_MG; 770 pc->ops->setfromoptions = PCSetFromOptions_MG; 771 pc->ops->view = PCView_MG; 772 773 pc->data = (void*)0; 774 PetscFunctionReturn(0); 775 } 776 EXTERN_C_END 777