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