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