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