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 #undef __FUNCT__ 10 #define __FUNCT__ "PCMGMCycle_Private" 11 PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 12 { 13 PC_MG *mg = (PC_MG*)pc->data; 14 PC_MG_Levels *mgc,*mglevels = *mglevelsin; 15 PetscErrorCode ierr; 16 PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 17 18 PetscFunctionBegin; 19 20 if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 21 ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 22 if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 23 if (mglevels->level) { /* not the coarsest grid */ 24 if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);} 25 ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 26 if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);} 27 28 /* if on finest level and have convergence criteria set */ 29 if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 30 PetscReal rnorm; 31 ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 32 if (rnorm <= mg->ttol) { 33 if (rnorm < mg->abstol) { 34 *reason = PCRICHARDSON_CONVERGED_ATOL; 35 ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr); 36 } else { 37 *reason = PCRICHARDSON_CONVERGED_RTOL; 38 ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr); 39 } 40 PetscFunctionReturn(0); 41 } 42 } 43 44 mgc = *(mglevelsin - 1); 45 if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 46 ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 47 if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 48 ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 49 while (cycles--) { 50 ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 51 } 52 if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 53 ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 54 if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 55 if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 56 ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 57 if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 58 } 59 PetscFunctionReturn(0); 60 } 61 62 #undef __FUNCT__ 63 #define __FUNCT__ "PCApplyRichardson_MG" 64 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason) 65 { 66 PC_MG *mg = (PC_MG*)pc->data; 67 PC_MG_Levels **mglevels = mg->levels; 68 PetscErrorCode ierr; 69 PetscInt levels = mglevels[0]->levels,i; 70 71 PetscFunctionBegin; 72 mglevels[levels-1]->b = b; 73 mglevels[levels-1]->x = x; 74 75 mg->rtol = rtol; 76 mg->abstol = abstol; 77 mg->dtol = dtol; 78 if (rtol) { 79 /* compute initial residual norm for relative convergence test */ 80 PetscReal rnorm; 81 if (zeroguess) { 82 ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 83 } else { 84 ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 85 ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 86 } 87 mg->ttol = PetscMax(rtol*rnorm,abstol); 88 } else if (abstol) { 89 mg->ttol = abstol; 90 } else { 91 mg->ttol = 0.0; 92 } 93 94 /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 95 stop prematurely do to small residual */ 96 for (i=1; i<levels; i++) { 97 ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 98 if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 99 ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 100 } 101 } 102 103 *reason = (PCRichardsonConvergedReason)0; 104 for (i=0; i<its; i++) { 105 ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 106 if (*reason) break; 107 } 108 if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 109 *outits = i; 110 PetscFunctionReturn(0); 111 } 112 113 #undef __FUNCT__ 114 #define __FUNCT__ "PCMGSetLevels" 115 /*@C 116 PCMGSetLevels - Sets the number of levels to use with MG. 117 Must be called before any other MG routine. 118 119 Collective on PC 120 121 Input Parameters: 122 + pc - the preconditioner context 123 . levels - the number of levels 124 - comms - optional communicators for each level; this is to allow solving the coarser problems 125 on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 126 127 Level: intermediate 128 129 Notes: 130 If the number of levels is one then the multigrid uses the -mg_levels prefix 131 for setting the level options rather than the -mg_coarse prefix. 132 133 .keywords: MG, set, levels, multigrid 134 135 .seealso: PCMGSetType(), PCMGGetLevels() 136 @*/ 137 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 138 { 139 PetscErrorCode ierr; 140 PC_MG *mg = (PC_MG*)pc->data; 141 MPI_Comm comm = ((PetscObject)pc)->comm; 142 PC_MG_Levels **mglevels; 143 PetscInt i; 144 PetscMPIInt size; 145 const char *prefix; 146 PC ipc; 147 148 PetscFunctionBegin; 149 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 150 if (mg->nlevels > -1) { 151 SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n make sure that you call PCMGSetLevels() before KSPSetFromOptions()"); 152 } 153 if (mg->levels) SETERRQ(PETSC_ERR_PLIB,"Internal error in PETSc, this array should not yet exist"); 154 155 mg->nlevels = levels; 156 mg->galerkin = PETSC_FALSE; 157 mg->galerkinused = PETSC_FALSE; 158 159 ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr); 160 ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 161 162 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 163 164 for (i=0; i<levels; i++) { 165 ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr); 166 mglevels[i]->level = i; 167 mglevels[i]->levels = levels; 168 mglevels[i]->cycles = PC_MG_CYCLE_V; 169 mg->default_smoothu = 1; 170 mg->default_smoothd = 1; 171 172 if (comms) comm = comms[i]; 173 ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 174 ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 175 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 176 ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 177 178 /* do special stuff for coarse grid */ 179 if (!i && levels > 1) { 180 ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 181 182 /* coarse solve is (redundant) LU by default */ 183 ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 184 ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 185 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 186 if (size > 1) { 187 ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 188 } else { 189 ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 190 } 191 192 } else { 193 char tprefix[128]; 194 sprintf(tprefix,"mg_levels_%d_",(int)i); 195 ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 196 } 197 ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr); 198 mglevels[i]->smoothu = mglevels[i]->smoothd; 199 mg->rtol = 0.0; 200 mg->abstol = 0.0; 201 mg->dtol = 0.0; 202 mg->ttol = 0.0; 203 mg->eventsmoothsetup = 0; 204 mg->eventsmoothsolve = 0; 205 mg->eventresidual = 0; 206 mg->eventinterprestrict = 0; 207 mg->cyclesperpcapply = 1; 208 } 209 mg->am = PC_MG_MULTIPLICATIVE; 210 mg->levels = mglevels; 211 pc->ops->applyrichardson = PCApplyRichardson_MG; 212 PetscFunctionReturn(0); 213 } 214 215 #undef __FUNCT__ 216 #define __FUNCT__ "PCDestroy_MG_Private" 217 PetscErrorCode PCDestroy_MG_Private(PC pc) 218 { 219 PC_MG *mg = (PC_MG*)pc->data; 220 PC_MG_Levels **mglevels = mg->levels; 221 PetscErrorCode ierr; 222 PetscInt i,n; 223 224 PetscFunctionBegin; 225 if (mglevels) { 226 n = mglevels[0]->levels; 227 for (i=0; i<n-1; i++) { 228 if (mglevels[i+1]->r) {ierr = VecDestroy(mglevels[i+1]->r);CHKERRQ(ierr);} 229 if (mglevels[i]->b) {ierr = VecDestroy(mglevels[i]->b);CHKERRQ(ierr);} 230 if (mglevels[i]->x) {ierr = VecDestroy(mglevels[i]->x);CHKERRQ(ierr);} 231 if (mglevels[i+1]->restrct) {ierr = MatDestroy(mglevels[i+1]->restrct);CHKERRQ(ierr);} 232 if (mglevels[i+1]->interpolate) {ierr = MatDestroy(mglevels[i+1]->interpolate);CHKERRQ(ierr);} 233 } 234 235 for (i=0; i<n; i++) { 236 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 237 ierr = KSPDestroy(mglevels[i]->smoothd);CHKERRQ(ierr); 238 } 239 ierr = KSPDestroy(mglevels[i]->smoothu);CHKERRQ(ierr); 240 ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 241 } 242 ierr = PetscFree(mglevels);CHKERRQ(ierr); 243 } 244 mg->nlevels = -1; 245 mg->levels = PETSC_NULL; 246 PetscFunctionReturn(0); 247 } 248 249 #undef __FUNCT__ 250 #define __FUNCT__ "PCDestroy_MG" 251 PetscErrorCode PCDestroy_MG(PC pc) 252 { 253 PC_MG *mg = (PC_MG*)pc->data; 254 PetscErrorCode ierr; 255 256 PetscFunctionBegin; 257 ierr = PCDestroy_MG_Private(pc);CHKERRQ(ierr); 258 ierr = PetscFree(mg);CHKERRQ(ierr); 259 PetscFunctionReturn(0); 260 } 261 262 263 264 EXTERN PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 265 EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 266 EXTERN PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 267 268 /* 269 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 270 or full cycle of multigrid. 271 272 Note: 273 A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 274 */ 275 #undef __FUNCT__ 276 #define __FUNCT__ "PCApply_MG" 277 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 278 { 279 PC_MG *mg = (PC_MG*)pc->data; 280 PC_MG_Levels **mglevels = mg->levels; 281 PetscErrorCode ierr; 282 PetscInt levels = mglevels[0]->levels,i; 283 284 PetscFunctionBegin; 285 mglevels[levels-1]->b = b; 286 mglevels[levels-1]->x = x; 287 if (mg->am == PC_MG_MULTIPLICATIVE) { 288 ierr = VecSet(x,0.0);CHKERRQ(ierr); 289 for (i=0; i<mg->cyclesperpcapply; i++) { 290 ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr); 291 } 292 } 293 else if (mg->am == PC_MG_ADDITIVE) { 294 ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 295 } 296 else if (mg->am == PC_MG_KASKADE) { 297 ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 298 } 299 else { 300 ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 301 } 302 PetscFunctionReturn(0); 303 } 304 305 306 #undef __FUNCT__ 307 #define __FUNCT__ "PCSetFromOptions_MG" 308 PetscErrorCode PCSetFromOptions_MG(PC pc) 309 { 310 PetscErrorCode ierr; 311 PetscInt m,levels = 1,cycles; 312 PetscTruth flg; 313 PC_MG *mg = (PC_MG*)pc->data; 314 PC_MG_Levels **mglevels = mg->levels; 315 PCMGType mgtype; 316 PCMGCycleType mgctype; 317 318 PetscFunctionBegin; 319 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 320 if (!mglevels) { 321 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 322 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 323 mglevels = mg->levels; 324 } 325 mgctype = (PCMGCycleType) mglevels[0]->cycles; 326 ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 327 if (flg) { 328 ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 329 }; 330 flg = PETSC_FALSE; 331 ierr = PetscOptionsTruth("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 332 if (flg) { 333 ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 334 } 335 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 336 if (flg) { 337 ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 338 } 339 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 340 if (flg) { 341 ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 342 } 343 mgtype = mg->am; 344 ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 345 if (flg) { 346 ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 347 } 348 if (mg->am == PC_MG_MULTIPLICATIVE) { 349 ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 350 if (flg) { 351 ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 352 } 353 } 354 flg = PETSC_FALSE; 355 ierr = PetscOptionsTruth("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 356 if (flg) { 357 PetscInt i; 358 char eventname[128]; 359 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 360 levels = mglevels[0]->levels; 361 for (i=0; i<levels; i++) { 362 sprintf(eventname,"MGSetup Level %d",(int)i); 363 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mg->eventsmoothsetup);CHKERRQ(ierr); 364 sprintf(eventname,"MGSmooth Level %d",(int)i); 365 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mg->eventsmoothsolve);CHKERRQ(ierr); 366 if (i) { 367 sprintf(eventname,"MGResid Level %d",(int)i); 368 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mg->eventresidual);CHKERRQ(ierr); 369 sprintf(eventname,"MGInterp Level %d",(int)i); 370 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mg->eventinterprestrict);CHKERRQ(ierr); 371 } 372 } 373 } 374 ierr = PetscOptionsTail();CHKERRQ(ierr); 375 PetscFunctionReturn(0); 376 } 377 378 const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 379 const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 380 381 #undef __FUNCT__ 382 #define __FUNCT__ "PCView_MG" 383 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 384 { 385 PC_MG *mg = (PC_MG*)pc->data; 386 PC_MG_Levels **mglevels = mg->levels; 387 PetscErrorCode ierr; 388 PetscInt levels = mglevels[0]->levels,i; 389 PetscTruth iascii; 390 391 PetscFunctionBegin; 392 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 393 if (iascii) { 394 ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,(mglevels[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr); 395 if (mg->am == PC_MG_MULTIPLICATIVE) { 396 ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 397 } 398 if (mg->galerkin) { 399 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 400 } 401 for (i=0; i<levels; i++) { 402 if (!i) { 403 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D presmooths=%D postsmooths=%D -----\n",i,mg->default_smoothd,mg->default_smoothu);CHKERRQ(ierr); 404 } else { 405 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr); 406 } 407 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 408 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 409 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 410 if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 411 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 412 } else if (i){ 413 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr); 414 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 415 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 416 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 417 } 418 } 419 } else { 420 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 421 } 422 PetscFunctionReturn(0); 423 } 424 425 /* 426 Calls setup for the KSP on each level 427 */ 428 #undef __FUNCT__ 429 #define __FUNCT__ "PCSetUp_MG" 430 PetscErrorCode PCSetUp_MG(PC pc) 431 { 432 PC_MG *mg = (PC_MG*)pc->data; 433 PC_MG_Levels **mglevels = mg->levels; 434 PetscErrorCode ierr; 435 PetscInt i,n = mglevels[0]->levels; 436 PC cpc,mpc; 437 PetscTruth preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset; 438 PetscViewerASCIIMonitor ascii; 439 PetscViewer viewer = PETSC_NULL; 440 MPI_Comm comm; 441 Mat dA,dB; 442 MatStructure uflag; 443 Vec tvec; 444 DM *dms; 445 446 PetscFunctionBegin; 447 448 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 449 /* so use those from global PC */ 450 /* Is this what we always want? What if user wants to keep old one? */ 451 ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 452 ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 453 ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr); 454 if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) { 455 ierr = PetscInfo(pc,"Using outer operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr); 456 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 457 } 458 459 if (pc->dm) { 460 ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 461 dms[n-1] = pc->dm; 462 for (i=n-2; i>-1; i--) { 463 ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr); 464 } 465 for (i=n-2; i>-1; i--) { 466 ierr = DMDestroy(dms[i]);CHKERRQ(ierr); 467 } 468 } 469 470 if (mg->galerkin) { 471 Mat B; 472 mg->galerkinused = PETSC_TRUE; 473 /* currently only handle case where mat and pmat are the same on coarser levels */ 474 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 475 if (!pc->setupcalled) { 476 for (i=n-2; i>-1; i--) { 477 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 478 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 479 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 480 dB = B; 481 } 482 ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr); 483 } else { 484 for (i=n-2; i>-1; i--) { 485 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 486 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 487 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 488 dB = B; 489 } 490 } 491 } 492 493 if (!pc->setupcalled) { 494 ierr = PetscOptionsGetTruth(0,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr); 495 496 for (i=0; i<n; i++) { 497 if (monitor) { 498 ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr); 499 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 500 ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 501 } 502 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 503 } 504 for (i=1; i<n; i++) { 505 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 506 if (monitor) { 507 ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr); 508 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 509 ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 510 } 511 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 512 } 513 } 514 for (i=1; i<n; i++) { 515 if (!mglevels[i]->residual) { 516 Mat mat; 517 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 518 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 519 } 520 if (mglevels[i]->restrct && !mglevels[i]->interpolate) { 521 ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr); 522 } 523 if (!mglevels[i]->restrct && mglevels[i]->interpolate) { 524 ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr); 525 } 526 #if defined(PETSC_USE_DEBUG) 527 if (!mglevels[i]->restrct || !mglevels[i]->interpolate) { 528 SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 529 } 530 #endif 531 } 532 for (i=0; i<n-1; i++) { 533 if (!mglevels[i]->b) { 534 Vec *vec; 535 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 536 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 537 ierr = VecDestroy(*vec);CHKERRQ(ierr); 538 ierr = PetscFree(vec);CHKERRQ(ierr); 539 } 540 if (!mglevels[i]->r && i) { 541 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 542 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 543 ierr = VecDestroy(tvec);CHKERRQ(ierr); 544 } 545 if (!mglevels[i]->x) { 546 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 547 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 548 ierr = VecDestroy(tvec);CHKERRQ(ierr); 549 } 550 } 551 if (n != 1 && !mglevels[n-1]->r) { 552 /* PCMGSetR() on the finest level if user did not supply it */ 553 Vec *vec; 554 ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 555 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 556 ierr = VecDestroy(*vec);CHKERRQ(ierr); 557 ierr = PetscFree(vec);CHKERRQ(ierr); 558 } 559 } 560 561 562 for (i=1; i<n; i++) { 563 if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 564 /* if doing only down then initial guess is zero */ 565 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 566 } 567 if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 568 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 569 if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 570 } 571 for (i=1; i<n; i++) { 572 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 573 Mat downmat,downpmat; 574 MatStructure matflag; 575 PetscTruth opsset; 576 577 /* check if operators have been set for up, if not use down operators to set them */ 578 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 579 if (!opsset) { 580 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 581 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 582 } 583 584 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 585 if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 586 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 587 if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 588 } 589 } 590 591 /* 592 If coarse solver is not direct method then DO NOT USE preonly 593 */ 594 ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 595 if (preonly) { 596 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 597 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 598 ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 599 if (!lu && !redundant && !cholesky) { 600 ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 601 } 602 } 603 604 if (!pc->setupcalled) { 605 if (monitor) { 606 ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr); 607 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 608 ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 609 } 610 ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 611 } 612 613 if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 614 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 615 if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 616 617 /* 618 Dump the interpolation/restriction matrices plus the 619 Jacobian/stiffness on each level. This allows Matlab users to 620 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 621 622 Only support one or the other at the same time. 623 */ 624 #if defined(PETSC_USE_SOCKET_VIEWER) 625 ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 626 if (dump) { 627 viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 628 } 629 dump = PETSC_FALSE; 630 #endif 631 ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 632 if (dump) { 633 viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 634 } 635 636 if (viewer) { 637 for (i=1; i<n; i++) { 638 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 639 } 640 for (i=0; i<n; i++) { 641 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 642 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 643 } 644 } 645 PetscFunctionReturn(0); 646 } 647 648 /* -------------------------------------------------------------------------------------*/ 649 650 #undef __FUNCT__ 651 #define __FUNCT__ "PCMGGetLevels" 652 /*@ 653 PCMGGetLevels - Gets the number of levels to use with MG. 654 655 Not Collective 656 657 Input Parameter: 658 . pc - the preconditioner context 659 660 Output parameter: 661 . levels - the number of levels 662 663 Level: advanced 664 665 .keywords: MG, get, levels, multigrid 666 667 .seealso: PCMGSetLevels() 668 @*/ 669 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels) 670 { 671 PC_MG *mg = (PC_MG*)pc->data; 672 673 PetscFunctionBegin; 674 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 675 PetscValidIntPointer(levels,2); 676 *levels = mg->nlevels; 677 PetscFunctionReturn(0); 678 } 679 680 #undef __FUNCT__ 681 #define __FUNCT__ "PCMGSetType" 682 /*@ 683 PCMGSetType - Determines the form of multigrid to use: 684 multiplicative, additive, full, or the Kaskade algorithm. 685 686 Collective on PC 687 688 Input Parameters: 689 + pc - the preconditioner context 690 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 691 PC_MG_FULL, PC_MG_KASKADE 692 693 Options Database Key: 694 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 695 additive, full, kaskade 696 697 Level: advanced 698 699 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 700 701 .seealso: PCMGSetLevels() 702 @*/ 703 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form) 704 { 705 PC_MG *mg = (PC_MG*)pc->data; 706 707 PetscFunctionBegin; 708 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 709 mg->am = form; 710 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 711 else pc->ops->applyrichardson = 0; 712 PetscFunctionReturn(0); 713 } 714 715 #undef __FUNCT__ 716 #define __FUNCT__ "PCMGSetCycleType" 717 /*@ 718 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 719 complicated cycling. 720 721 Collective on PC 722 723 Input Parameters: 724 + pc - the multigrid context 725 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 726 727 Options Database Key: 728 $ -pc_mg_cycle_type v or w 729 730 Level: advanced 731 732 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 733 734 .seealso: PCMGSetCycleTypeOnLevel() 735 @*/ 736 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n) 737 { 738 PC_MG *mg = (PC_MG*)pc->data; 739 PC_MG_Levels **mglevels = mg->levels; 740 PetscInt i,levels; 741 742 PetscFunctionBegin; 743 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 744 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 745 levels = mglevels[0]->levels; 746 747 for (i=0; i<levels; i++) { 748 mglevels[i]->cycles = n; 749 } 750 PetscFunctionReturn(0); 751 } 752 753 #undef __FUNCT__ 754 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 755 /*@ 756 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 757 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 758 759 Collective on PC 760 761 Input Parameters: 762 + pc - the multigrid context 763 - n - number of cycles (default is 1) 764 765 Options Database Key: 766 $ -pc_mg_multiplicative_cycles n 767 768 Level: advanced 769 770 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 771 772 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 773 774 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 775 @*/ 776 PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 777 { 778 PC_MG *mg = (PC_MG*)pc->data; 779 PC_MG_Levels **mglevels = mg->levels; 780 PetscInt i,levels; 781 782 PetscFunctionBegin; 783 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 784 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 785 levels = mglevels[0]->levels; 786 787 for (i=0; i<levels; i++) { 788 mg->cyclesperpcapply = n; 789 } 790 PetscFunctionReturn(0); 791 } 792 793 #undef __FUNCT__ 794 #define __FUNCT__ "PCMGSetGalerkin" 795 /*@ 796 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 797 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 798 799 Collective on PC 800 801 Input Parameters: 802 . pc - the multigrid context 803 804 Options Database Key: 805 $ -pc_mg_galerkin 806 807 Level: intermediate 808 809 .keywords: MG, set, Galerkin 810 811 .seealso: PCMGGetGalerkin() 812 813 @*/ 814 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 815 { 816 PC_MG *mg = (PC_MG*)pc->data; 817 818 PetscFunctionBegin; 819 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 820 mg->galerkin = PETSC_TRUE; 821 PetscFunctionReturn(0); 822 } 823 824 #undef __FUNCT__ 825 #define __FUNCT__ "PCMGGetGalerkin" 826 /*@ 827 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 828 A_i-1 = r_i * A_i * r_i^t 829 830 Not Collective 831 832 Input Parameter: 833 . pc - the multigrid context 834 835 Output Parameter: 836 . gelerkin - PETSC_TRUE or PETSC_FALSE 837 838 Options Database Key: 839 $ -pc_mg_galerkin 840 841 Level: intermediate 842 843 .keywords: MG, set, Galerkin 844 845 .seealso: PCMGSetGalerkin() 846 847 @*/ 848 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 849 { 850 PC_MG *mg = (PC_MG*)pc->data; 851 852 PetscFunctionBegin; 853 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 854 *galerkin = mg->galerkin; 855 PetscFunctionReturn(0); 856 } 857 858 #undef __FUNCT__ 859 #define __FUNCT__ "PCMGSetNumberSmoothDown" 860 /*@ 861 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 862 use on all levels. Use PCMGGetSmootherDown() to set different 863 pre-smoothing steps on different levels. 864 865 Collective on PC 866 867 Input Parameters: 868 + mg - the multigrid context 869 - n - the number of smoothing steps 870 871 Options Database Key: 872 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 873 874 Level: advanced 875 876 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 877 878 .seealso: PCMGSetNumberSmoothUp() 879 @*/ 880 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n) 881 { 882 PC_MG *mg = (PC_MG*)pc->data; 883 PC_MG_Levels **mglevels = mg->levels; 884 PetscErrorCode ierr; 885 PetscInt i,levels; 886 887 PetscFunctionBegin; 888 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 889 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 890 levels = mglevels[0]->levels; 891 892 for (i=1; i<levels; i++) { 893 /* make sure smoother up and down are different */ 894 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 895 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 896 mg->default_smoothd = n; 897 } 898 PetscFunctionReturn(0); 899 } 900 901 #undef __FUNCT__ 902 #define __FUNCT__ "PCMGSetNumberSmoothUp" 903 /*@ 904 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 905 on all levels. Use PCMGGetSmootherUp() to set different numbers of 906 post-smoothing steps on different levels. 907 908 Collective on PC 909 910 Input Parameters: 911 + mg - the multigrid context 912 - n - the number of smoothing steps 913 914 Options Database Key: 915 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 916 917 Level: advanced 918 919 Note: this does not set a value on the coarsest grid, since we assume that 920 there is no separate smooth up on the coarsest grid. 921 922 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 923 924 .seealso: PCMGSetNumberSmoothDown() 925 @*/ 926 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 927 { 928 PC_MG *mg = (PC_MG*)pc->data; 929 PC_MG_Levels **mglevels = mg->levels; 930 PetscErrorCode ierr; 931 PetscInt i,levels; 932 933 PetscFunctionBegin; 934 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 935 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 936 levels = mglevels[0]->levels; 937 938 for (i=1; i<levels; i++) { 939 /* make sure smoother up and down are different */ 940 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 941 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 942 mg->default_smoothu = n; 943 } 944 PetscFunctionReturn(0); 945 } 946 947 /* ----------------------------------------------------------------------------------------*/ 948 949 /*MC 950 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 951 information about the coarser grid matrices and restriction/interpolation operators. 952 953 Options Database Keys: 954 + -pc_mg_levels <nlevels> - number of levels including finest 955 . -pc_mg_cycles v or w 956 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 957 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 958 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 959 . -pc_mg_log - log information about time spent on each level of the solver 960 . -pc_mg_monitor - print information on the multigrid convergence 961 . -pc_mg_galerkin - use Galerkin process to compute coarser operators 962 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 963 to the Socket viewer for reading from Matlab. 964 965 Notes: 966 967 Level: intermediate 968 969 Concepts: multigrid/multilevel 970 971 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 972 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 973 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 974 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 975 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 976 M*/ 977 978 EXTERN_C_BEGIN 979 #undef __FUNCT__ 980 #define __FUNCT__ "PCCreate_MG" 981 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 982 { 983 PC_MG *mg; 984 PetscErrorCode ierr; 985 986 PetscFunctionBegin; 987 ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 988 pc->data = (void*)mg; 989 mg->nlevels = -1; 990 991 pc->ops->apply = PCApply_MG; 992 pc->ops->setup = PCSetUp_MG; 993 pc->ops->destroy = PCDestroy_MG; 994 pc->ops->setfromoptions = PCSetFromOptions_MG; 995 pc->ops->view = PCView_MG; 996 PetscFunctionReturn(0); 997 } 998 EXTERN_C_END 999