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