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 int 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 int MGACycle_Private(MG*); 148 EXTERN int MGFCycle_Private(MG*); 149 EXTERN int 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 int ierr,levels = mg[0]->levels; 191 PetscTruth converged = PETSC_FALSE; 192 193 PetscFunctionBegin; 194 mg[levels-1]->b = b; 195 mg[levels-1]->x = x; 196 197 mg[levels-1]->rtol = rtol; 198 mg[levels-1]->atol = atol; 199 mg[levels-1]->dtol = dtol; 200 if (rtol) { 201 /* compute initial residual norm for relative convergence test */ 202 PetscReal rnorm; 203 ierr = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr); 204 ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 205 mg[levels-1]->ttol = PetscMax(rtol*rnorm,atol); 206 } else if (atol) { 207 mg[levels-1]->ttol = atol; 208 } else { 209 mg[levels-1]->ttol = 0.0; 210 } 211 212 while (its-- && !converged) { 213 ierr = MGMCycle_Private(mg+levels-1,&converged);CHKERRQ(ierr); 214 } 215 PetscFunctionReturn(0); 216 } 217 218 #undef __FUNCT__ 219 #define __FUNCT__ "PCSetFromOptions_MG" 220 static int PCSetFromOptions_MG(PC pc) 221 { 222 int ierr,indx,m,levels = 1; 223 PetscTruth flg; 224 const char *type[] = {"additive","multiplicative","full","cascade","kascade"}; 225 226 PetscFunctionBegin; 227 228 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 229 if (!pc->data) { 230 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","MGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 231 ierr = MGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 232 } 233 ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 234 if (flg) { 235 ierr = MGSetCycles(pc,m);CHKERRQ(ierr); 236 } 237 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 238 if (flg) { 239 ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 240 } 241 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 242 if (flg) { 243 ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 244 } 245 ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr); 246 if (flg) { 247 MGType mg = (MGType) 0; 248 switch (indx) { 249 case 0: 250 mg = MGADDITIVE; 251 break; 252 case 1: 253 mg = MGMULTIPLICATIVE; 254 break; 255 case 2: 256 mg = MGFULL; 257 break; 258 case 3: 259 mg = MGKASKADE; 260 break; 261 case 4: 262 mg = MGKASKADE; 263 break; 264 } 265 ierr = MGSetType(pc,mg);CHKERRQ(ierr); 266 } 267 ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr); 268 if (flg) { 269 MG *mg = (MG*)pc->data; 270 int i; 271 char eventname[128]; 272 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 273 levels = mg[0]->levels; 274 for (i=0; i<levels; i++) { 275 sprintf(eventname,"MSetup Level %d",i); 276 ierr = PetscLogEventRegister(&mg[i]->eventsetup,eventname,pc->cookie);CHKERRQ(ierr); 277 sprintf(eventname,"MGSolve Level %d",i); 278 ierr = PetscLogEventRegister(&mg[i]->eventsolve,eventname,pc->cookie);CHKERRQ(ierr); 279 } 280 } 281 ierr = PetscOptionsTail();CHKERRQ(ierr); 282 PetscFunctionReturn(0); 283 } 284 285 #undef __FUNCT__ 286 #define __FUNCT__ "PCView_MG" 287 static int PCView_MG(PC pc,PetscViewer viewer) 288 { 289 MG *mg = (MG*)pc->data; 290 int ierr,levels = mg[0]->levels,i; 291 const char *cstring; 292 PetscTruth iascii; 293 294 PetscFunctionBegin; 295 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 296 if (iascii) { 297 if (mg[0]->am == MGMULTIPLICATIVE) cstring = "multiplicative"; 298 else if (mg[0]->am == MGADDITIVE) cstring = "additive"; 299 else if (mg[0]->am == MGFULL) cstring = "full"; 300 else if (mg[0]->am == MGKASKADE) cstring = "Kaskade"; 301 else cstring = "unknown"; 302 ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%d cycles=%d, pre-smooths=%d, post-smooths=%d\n", 303 cstring,levels,mg[0]->cycles,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr); 304 for (i=0; i<levels; i++) { 305 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %d -------------------------------\n",i);CHKERRQ(ierr); 306 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 307 ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr); 308 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 309 if (mg[i]->smoothd == mg[i]->smoothu) { 310 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 311 } else { 312 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %d -------------------------------\n",i);CHKERRQ(ierr); 313 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 314 ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr); 315 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 316 } 317 } 318 } else { 319 SETERRQ1(1,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 320 } 321 PetscFunctionReturn(0); 322 } 323 324 /* 325 Calls setup for the KSP on each level 326 */ 327 #undef __FUNCT__ 328 #define __FUNCT__ "PCSetUp_MG" 329 static int PCSetUp_MG(PC pc) 330 { 331 MG *mg = (MG*)pc->data; 332 int ierr,i,n = mg[0]->levels; 333 PC cpc; 334 PetscTruth preonly,lu,redundant,monitor = PETSC_FALSE,dump; 335 PetscViewer ascii; 336 MPI_Comm comm; 337 338 PetscFunctionBegin; 339 if (!pc->setupcalled) { 340 ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr); 341 342 for (i=1; i<n; i++) { 343 if (mg[i]->smoothd) { 344 if (monitor) { 345 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr); 346 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 347 ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr); 348 ierr = KSPSetMonitor(mg[i]->smoothd,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 349 } 350 ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr); 351 } 352 } 353 for (i=0; i<n; i++) { 354 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 355 if (monitor) { 356 ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr); 357 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 358 ierr = PetscViewerASCIISetTab(ascii,n-i);CHKERRQ(ierr); 359 ierr = KSPSetMonitor(mg[i]->smoothu,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 360 } 361 ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr); 362 } 363 } 364 } 365 366 for (i=1; i<n; i++) { 367 if (mg[i]->smoothd) { 368 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 369 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 370 ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 371 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 372 } 373 } 374 for (i=0; i<n; i++) { 375 if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 376 ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 377 if (mg[i]->eventsetup) {ierr = PetscLogEventBegin(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 378 ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 379 if (mg[i]->eventsetup) {ierr = PetscLogEventEnd(mg[i]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 380 } 381 } 382 383 /* 384 If coarse solver is not direct method then DO NOT USE preonly 385 */ 386 ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 387 if (preonly) { 388 ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 389 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 390 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 391 if (!lu && !redundant) { 392 ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 393 } 394 } 395 396 if (!pc->setupcalled) { 397 if (monitor) { 398 ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 399 ierr = PetscViewerASCIIOpen(comm,"stdout",&ascii);CHKERRQ(ierr); 400 ierr = PetscViewerASCIISetTab(ascii,n);CHKERRQ(ierr); 401 ierr = KSPSetMonitor(mg[0]->smoothd,KSPDefaultMonitor,ascii,(int(*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 402 } 403 ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr); 404 } 405 406 if (mg[0]->eventsetup) {ierr = PetscLogEventBegin(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 407 ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr); 408 if (mg[0]->eventsetup) {ierr = PetscLogEventEnd(mg[0]->eventsetup,0,0,0,0);CHKERRQ(ierr);} 409 410 /* 411 Dump the interpolation/restriction matrices to matlab plus the 412 Jacobian/stiffness on each level. This allows Matlab users to 413 easily check if the Galerkin condition A_c = R A_f R^T is satisfied */ 414 ierr = PetscOptionsHasName(pc->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr); 415 if (dump) { 416 for (i=1; i<n; i++) { 417 ierr = MatView(mg[i]->restrct,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr); 418 } 419 for (i=0; i<n; i++) { 420 ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 421 ierr = MatView(pc->mat,PETSC_VIEWER_SOCKET_(pc->comm));CHKERRQ(ierr); 422 } 423 } 424 425 PetscFunctionReturn(0); 426 } 427 428 /* -------------------------------------------------------------------------------------*/ 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "MGSetLevels" 432 /*@C 433 MGSetLevels - Sets the number of levels to use with MG. 434 Must be called before any other MG routine. 435 436 Collective on PC 437 438 Input Parameters: 439 + pc - the preconditioner context 440 . levels - the number of levels 441 - comms - optional communicators for each level; this is to allow solving the coarser problems 442 on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 443 444 Level: intermediate 445 446 Notes: 447 If the number of levels is one then the multigrid uses the -mg_levels prefix 448 for setting the level options rather than the -mg_coarse prefix. 449 450 .keywords: MG, set, levels, multigrid 451 452 .seealso: MGSetType(), MGGetLevels() 453 @*/ 454 int MGSetLevels(PC pc,int levels,MPI_Comm *comms) 455 { 456 int ierr; 457 MG *mg; 458 459 PetscFunctionBegin; 460 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 461 462 if (pc->data) { 463 SETERRQ(1,"Number levels already set for MG\n\ 464 make sure that you call MGSetLevels() before KSPSetFromOptions()"); 465 } 466 ierr = MGCreate_Private(pc->comm,levels,pc,comms,&mg);CHKERRQ(ierr); 467 mg[0]->am = MGMULTIPLICATIVE; 468 pc->data = (void*)mg; 469 pc->ops->applyrichardson = PCApplyRichardson_MG; 470 PetscFunctionReturn(0); 471 } 472 473 #undef __FUNCT__ 474 #define __FUNCT__ "MGGetLevels" 475 /*@ 476 MGGetLevels - Gets the number of levels to use with MG. 477 478 Not Collective 479 480 Input Parameter: 481 . pc - the preconditioner context 482 483 Output parameter: 484 . levels - the number of levels 485 486 Level: advanced 487 488 .keywords: MG, get, levels, multigrid 489 490 .seealso: MGSetLevels() 491 @*/ 492 int MGGetLevels(PC pc,int *levels) 493 { 494 MG *mg; 495 496 PetscFunctionBegin; 497 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 498 PetscValidIntPointer(levels,2); 499 500 mg = (MG*)pc->data; 501 *levels = mg[0]->levels; 502 PetscFunctionReturn(0); 503 } 504 505 #undef __FUNCT__ 506 #define __FUNCT__ "MGSetType" 507 /*@ 508 MGSetType - Determines the form of multigrid to use: 509 multiplicative, additive, full, or the Kaskade algorithm. 510 511 Collective on PC 512 513 Input Parameters: 514 + pc - the preconditioner context 515 - form - multigrid form, one of MGMULTIPLICATIVE, MGADDITIVE, 516 MGFULL, MGKASKADE 517 518 Options Database Key: 519 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 520 additive, full, kaskade 521 522 Level: advanced 523 524 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 525 526 .seealso: MGSetLevels() 527 @*/ 528 int MGSetType(PC pc,MGType form) 529 { 530 MG *mg; 531 532 PetscFunctionBegin; 533 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 534 mg = (MG*)pc->data; 535 536 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 537 mg[0]->am = form; 538 if (form == MGMULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 539 else pc->ops->applyrichardson = 0; 540 PetscFunctionReturn(0); 541 } 542 543 #undef __FUNCT__ 544 #define __FUNCT__ "MGSetCycles" 545 /*@ 546 MGSetCycles - Sets the type cycles to use. Use MGSetCyclesOnLevel() for more 547 complicated cycling. 548 549 Collective on PC 550 551 Input Parameters: 552 + mg - the multigrid context 553 - n - the number of cycles 554 555 Options Database Key: 556 $ -pc_mg_cycles n - 1 denotes a V-cycle; 2 denotes a W-cycle. 557 558 Level: advanced 559 560 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 561 562 .seealso: MGSetCyclesOnLevel() 563 @*/ 564 int MGSetCycles(PC pc,int n) 565 { 566 MG *mg; 567 int i,levels; 568 569 PetscFunctionBegin; 570 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 571 mg = (MG*)pc->data; 572 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 573 levels = mg[0]->levels; 574 575 for (i=0; i<levels; i++) { 576 mg[i]->cycles = n; 577 } 578 PetscFunctionReturn(0); 579 } 580 581 #undef __FUNCT__ 582 #define __FUNCT__ "MGCheck" 583 /*@ 584 MGCheck - Checks that all components of the MG structure have 585 been set. 586 587 Collective on PC 588 589 Input Parameters: 590 . mg - the MG structure 591 592 Level: advanced 593 594 .keywords: MG, check, set, multigrid 595 @*/ 596 int MGCheck(PC pc) 597 { 598 MG *mg; 599 int i,n,count = 0; 600 601 PetscFunctionBegin; 602 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 603 mg = (MG*)pc->data; 604 605 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 606 607 n = mg[0]->levels; 608 609 for (i=1; i<n; i++) { 610 if (!mg[i]->restrct) { 611 (*PetscErrorPrintf)("No restrict set level %d \n",n-i); count++; 612 } 613 if (!mg[i]->interpolate) { 614 (*PetscErrorPrintf)("No interpolate set level %d \n",n-i); count++; 615 } 616 if (!mg[i]->residual) { 617 (*PetscErrorPrintf)("No residual set level %d \n",n-i); count++; 618 } 619 if (!mg[i]->smoothu) { 620 (*PetscErrorPrintf)("No smoothup set level %d \n",n-i); count++; 621 } 622 if (!mg[i]->smoothd) { 623 (*PetscErrorPrintf)("No smoothdown set level %d \n",n-i); count++; 624 } 625 if (!mg[i]->r) { 626 (*PetscErrorPrintf)("No r set level %d \n",n-i); count++; 627 } 628 if (!mg[i-1]->x) { 629 (*PetscErrorPrintf)("No x set level %d \n",n-i); count++; 630 } 631 if (!mg[i-1]->b) { 632 (*PetscErrorPrintf)("No b set level %d \n",n-i); count++; 633 } 634 } 635 PetscFunctionReturn(count); 636 } 637 638 639 #undef __FUNCT__ 640 #define __FUNCT__ "MGSetNumberSmoothDown" 641 /*@ 642 MGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 643 use on all levels. Use MGGetSmootherDown() to set different 644 pre-smoothing steps on different levels. 645 646 Collective on PC 647 648 Input Parameters: 649 + mg - the multigrid context 650 - n - the number of smoothing steps 651 652 Options Database Key: 653 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 654 655 Level: advanced 656 657 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 658 659 .seealso: MGSetNumberSmoothUp() 660 @*/ 661 int MGSetNumberSmoothDown(PC pc,int n) 662 { 663 MG *mg; 664 int i,levels,ierr; 665 666 PetscFunctionBegin; 667 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 668 mg = (MG*)pc->data; 669 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 670 levels = mg[0]->levels; 671 672 for (i=0; i<levels; i++) { 673 /* make sure smoother up and down are different */ 674 ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 675 ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 676 mg[i]->default_smoothd = n; 677 } 678 PetscFunctionReturn(0); 679 } 680 681 #undef __FUNCT__ 682 #define __FUNCT__ "MGSetNumberSmoothUp" 683 /*@ 684 MGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 685 on all levels. Use MGGetSmootherUp() to set different numbers of 686 post-smoothing steps on different levels. 687 688 Collective on PC 689 690 Input Parameters: 691 + mg - the multigrid context 692 - n - the number of smoothing steps 693 694 Options Database Key: 695 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 696 697 Level: advanced 698 699 Note: this does not set a value on the coarsest grid, since we assume that 700 there is no seperate smooth up on the coarsest grid. If you want to have a 701 seperate smooth up on the coarsest grid then call MGGetSmoothUp(pc,0,&ksp); 702 703 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 704 705 .seealso: MGSetNumberSmoothDown() 706 @*/ 707 int MGSetNumberSmoothUp(PC pc,int n) 708 { 709 MG *mg; 710 int i,levels,ierr; 711 712 PetscFunctionBegin; 713 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 714 mg = (MG*)pc->data; 715 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 716 levels = mg[0]->levels; 717 718 for (i=1; i<levels; i++) { 719 /* make sure smoother up and down are different */ 720 ierr = MGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 721 ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 722 mg[i]->default_smoothu = n; 723 } 724 PetscFunctionReturn(0); 725 } 726 727 /* ----------------------------------------------------------------------------------------*/ 728 729 /*MC 730 PCMG - Use geometric multigrid preconditioning. This preconditioner requires you provide additional 731 information about the coarser grid matrices and restriction/interpolation operators. 732 733 Options Database Keys: 734 + -pc_mg_levels <nlevels> - number of levels including finest 735 . -pc_mg_cycles 1 or 2 - for V or W-cycle 736 . -pc_mg_smoothup <n> - number of smoothing steps before interpolation 737 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 738 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 739 . -pc_mg_log - log information about time spent on each level of the solver 740 . -pc_mg_monitor - print information on the multigrid convergence 741 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 742 to the Socket viewer for reading from Matlab. 743 744 Notes: 745 746 Level: intermediate 747 748 Concepts: multigrid 749 750 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 751 MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(), 752 MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(), 753 MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(), 754 MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR() 755 M*/ 756 757 EXTERN_C_BEGIN 758 #undef __FUNCT__ 759 #define __FUNCT__ "PCCreate_MG" 760 int PCCreate_MG(PC pc) 761 { 762 PetscFunctionBegin; 763 pc->ops->apply = PCApply_MG; 764 pc->ops->setup = PCSetUp_MG; 765 pc->ops->destroy = PCDestroy_MG; 766 pc->ops->setfromoptions = PCSetFromOptions_MG; 767 pc->ops->view = PCView_MG; 768 769 pc->data = (void*)0; 770 PetscFunctionReturn(0); 771 } 772 EXTERN_C_END 773