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