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