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_COOKIE,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 157 ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr); 158 ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 159 160 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 161 162 for (i=0; i<levels; i++) { 163 ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr); 164 mglevels[i]->level = i; 165 mglevels[i]->levels = levels; 166 mglevels[i]->cycles = PC_MG_CYCLE_V; 167 mglevels[i]->galerkin = PETSC_FALSE; 168 mglevels[i]->galerkinused = PETSC_FALSE; 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" 217 PetscErrorCode PCDestroy_MG(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 ierr = PetscFree(mg);CHKERRQ(ierr); 245 PetscFunctionReturn(0); 246 } 247 248 249 250 EXTERN PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 251 EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 252 EXTERN PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 253 254 /* 255 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 256 or full cycle of multigrid. 257 258 Note: 259 A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 260 */ 261 #undef __FUNCT__ 262 #define __FUNCT__ "PCApply_MG" 263 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 264 { 265 PC_MG *mg = (PC_MG*)pc->data; 266 PC_MG_Levels **mglevels = mg->levels; 267 PetscErrorCode ierr; 268 PetscInt levels = mglevels[0]->levels,i; 269 270 PetscFunctionBegin; 271 mglevels[levels-1]->b = b; 272 mglevels[levels-1]->x = x; 273 if (mg->am == PC_MG_MULTIPLICATIVE) { 274 ierr = VecSet(x,0.0);CHKERRQ(ierr); 275 for (i=0; i<mg->cyclesperpcapply; i++) { 276 ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr); 277 } 278 } 279 else if (mg->am == PC_MG_ADDITIVE) { 280 ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 281 } 282 else if (mg->am == PC_MG_KASKADE) { 283 ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 284 } 285 else { 286 ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 287 } 288 PetscFunctionReturn(0); 289 } 290 291 292 #undef __FUNCT__ 293 #define __FUNCT__ "PCSetFromOptions_MG" 294 PetscErrorCode PCSetFromOptions_MG(PC pc) 295 { 296 PetscErrorCode ierr; 297 PetscInt m,levels = 1,cycles; 298 PetscTruth flg; 299 PC_MG *mg = (PC_MG*)pc->data; 300 PC_MG_Levels **mglevels = mg->levels; 301 PCMGType mgtype; 302 PCMGCycleType mgctype; 303 304 PetscFunctionBegin; 305 ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 306 if (!pc->data) { 307 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 308 ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 309 mglevels = mg->levels; 310 } 311 mgctype = (PCMGCycleType) mglevels[0]->cycles; 312 ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 313 if (flg) { 314 ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 315 }; 316 flg = PETSC_FALSE; 317 ierr = PetscOptionsTruth("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 318 if (flg) { 319 ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 320 } 321 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 322 if (flg) { 323 ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 324 } 325 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 326 if (flg) { 327 ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 328 } 329 mgtype = mg->am; 330 ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 331 if (flg) { 332 ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 333 } 334 if (mg->am == PC_MG_MULTIPLICATIVE) { 335 ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 336 if (flg) { 337 ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 338 } 339 } 340 flg = PETSC_FALSE; 341 ierr = PetscOptionsTruth("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 342 if (flg) { 343 PetscInt i; 344 char eventname[128]; 345 if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 346 levels = mglevels[0]->levels; 347 for (i=0; i<levels; i++) { 348 sprintf(eventname,"MGSetup Level %d",(int)i); 349 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg->eventsmoothsetup);CHKERRQ(ierr); 350 sprintf(eventname,"MGSmooth Level %d",(int)i); 351 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg->eventsmoothsolve);CHKERRQ(ierr); 352 if (i) { 353 sprintf(eventname,"MGResid Level %d",(int)i); 354 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg->eventresidual);CHKERRQ(ierr); 355 sprintf(eventname,"MGInterp Level %d",(int)i); 356 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg->eventinterprestrict);CHKERRQ(ierr); 357 } 358 } 359 } 360 ierr = PetscOptionsTail();CHKERRQ(ierr); 361 PetscFunctionReturn(0); 362 } 363 364 const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 365 const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 366 367 #undef __FUNCT__ 368 #define __FUNCT__ "PCView_MG" 369 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 370 { 371 PC_MG *mg = (PC_MG*)pc->data; 372 PC_MG_Levels **mglevels = mg->levels; 373 PetscErrorCode ierr; 374 PetscInt levels = mglevels[0]->levels,i; 375 PetscTruth iascii; 376 377 PetscFunctionBegin; 378 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 379 if (iascii) { 380 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); 381 if (mg->am == PC_MG_MULTIPLICATIVE) { 382 ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 383 } 384 if (mglevels[0]->galerkin) { 385 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 386 } 387 for (i=0; i<levels; i++) { 388 if (!i) { 389 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D presmooths=%D postsmooths=%D -----\n",i,mg->default_smoothd,mg->default_smoothu);CHKERRQ(ierr); 390 } else { 391 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr); 392 } 393 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 394 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 395 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 396 if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 397 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 398 } else if (i){ 399 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr); 400 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 401 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 402 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 403 } 404 } 405 } else { 406 SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 407 } 408 PetscFunctionReturn(0); 409 } 410 411 /* 412 Calls setup for the KSP on each level 413 */ 414 #undef __FUNCT__ 415 #define __FUNCT__ "PCSetUp_MG" 416 PetscErrorCode PCSetUp_MG(PC pc) 417 { 418 PC_MG *mg = (PC_MG*)pc->data; 419 PC_MG_Levels **mglevels = mg->levels; 420 PetscErrorCode ierr; 421 PetscInt i,n = mglevels[0]->levels; 422 PC cpc,mpc; 423 PetscTruth preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset; 424 PetscViewerASCIIMonitor ascii; 425 PetscViewer viewer = PETSC_NULL; 426 MPI_Comm comm; 427 Mat dA,dB; 428 MatStructure uflag; 429 Vec tvec; 430 431 PetscFunctionBegin; 432 433 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 434 /* so use those from global PC */ 435 /* Is this what we always want? What if user wants to keep old one? */ 436 ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 437 ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 438 ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr); 439 if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) { 440 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); 441 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 442 } 443 444 if (mglevels[0]->galerkin) { 445 Mat B; 446 mglevels[0]->galerkinused = PETSC_TRUE; 447 /* currently only handle case where mat and pmat are the same on coarser levels */ 448 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 449 if (!pc->setupcalled) { 450 for (i=n-2; i>-1; i--) { 451 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 452 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 453 if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 454 dB = B; 455 } 456 ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr); 457 } else { 458 for (i=n-2; i>-1; i--) { 459 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 460 ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 461 ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 462 dB = B; 463 } 464 } 465 } 466 467 if (!pc->setupcalled) { 468 ierr = PetscOptionsGetTruth(0,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr); 469 470 for (i=0; i<n; i++) { 471 if (monitor) { 472 ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr); 473 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 474 ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 475 } 476 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 477 } 478 for (i=1; i<n; i++) { 479 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 480 if (monitor) { 481 ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr); 482 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 483 ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 484 } 485 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 486 } 487 } 488 for (i=1; i<n; i++) { 489 if (!mglevels[i]->residual) { 490 Mat mat; 491 ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 492 ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 493 } 494 if (mglevels[i]->restrct && !mglevels[i]->interpolate) { 495 ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr); 496 } 497 if (!mglevels[i]->restrct && mglevels[i]->interpolate) { 498 ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr); 499 } 500 #if defined(PETSC_USE_DEBUG) 501 if (!mglevels[i]->restrct || !mglevels[i]->interpolate) { 502 SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 503 } 504 #endif 505 } 506 for (i=0; i<n-1; i++) { 507 if (!mglevels[i]->b) { 508 Vec *vec; 509 ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 510 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 511 ierr = VecDestroy(*vec);CHKERRQ(ierr); 512 ierr = PetscFree(vec);CHKERRQ(ierr); 513 } 514 if (!mglevels[i]->r && i) { 515 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 516 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 517 ierr = VecDestroy(tvec);CHKERRQ(ierr); 518 } 519 if (!mglevels[i]->x) { 520 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 521 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 522 ierr = VecDestroy(tvec);CHKERRQ(ierr); 523 } 524 } 525 if (n != 1 && !mglevels[n-1]->r) { 526 /* PCMGSetR() on the finest level if user did not supply it */ 527 Vec *vec; 528 ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 529 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 530 ierr = VecDestroy(*vec);CHKERRQ(ierr); 531 ierr = PetscFree(vec);CHKERRQ(ierr); 532 } 533 } 534 535 536 for (i=1; i<n; i++) { 537 if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 538 /* if doing only down then initial guess is zero */ 539 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 540 } 541 if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 542 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 543 if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 544 } 545 for (i=1; i<n; i++) { 546 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 547 Mat downmat,downpmat; 548 MatStructure matflag; 549 PetscTruth opsset; 550 551 /* check if operators have been set for up, if not use down operators to set them */ 552 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 553 if (!opsset) { 554 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 555 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 556 } 557 558 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 559 if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 560 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 561 if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 562 } 563 } 564 565 /* 566 If coarse solver is not direct method then DO NOT USE preonly 567 */ 568 ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 569 if (preonly) { 570 ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 571 ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 572 ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 573 if (!lu && !redundant && !cholesky) { 574 ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 575 } 576 } 577 578 if (!pc->setupcalled) { 579 if (monitor) { 580 ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr); 581 ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 582 ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 583 } 584 ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 585 } 586 587 if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 588 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 589 if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 590 591 /* 592 Dump the interpolation/restriction matrices plus the 593 Jacobian/stiffness on each level. This allows Matlab users to 594 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 595 596 Only support one or the other at the same time. 597 */ 598 #if defined(PETSC_USE_SOCKET_VIEWER) 599 ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 600 if (dump) { 601 viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 602 } 603 dump = PETSC_FALSE; 604 #endif 605 ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 606 if (dump) { 607 viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 608 } 609 610 if (viewer) { 611 for (i=1; i<n; i++) { 612 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 613 } 614 for (i=0; i<n; i++) { 615 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 616 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 617 } 618 } 619 PetscFunctionReturn(0); 620 } 621 622 /* -------------------------------------------------------------------------------------*/ 623 624 #undef __FUNCT__ 625 #define __FUNCT__ "PCMGGetLevels" 626 /*@ 627 PCMGGetLevels - Gets the number of levels to use with MG. 628 629 Not Collective 630 631 Input Parameter: 632 . pc - the preconditioner context 633 634 Output parameter: 635 . levels - the number of levels 636 637 Level: advanced 638 639 .keywords: MG, get, levels, multigrid 640 641 .seealso: PCMGSetLevels() 642 @*/ 643 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels) 644 { 645 PC_MG *mg = (PC_MG*)pc->data; 646 647 PetscFunctionBegin; 648 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 649 PetscValidIntPointer(levels,2); 650 *levels = mg->nlevels; 651 PetscFunctionReturn(0); 652 } 653 654 #undef __FUNCT__ 655 #define __FUNCT__ "PCMGSetType" 656 /*@ 657 PCMGSetType - Determines the form of multigrid to use: 658 multiplicative, additive, full, or the Kaskade algorithm. 659 660 Collective on PC 661 662 Input Parameters: 663 + pc - the preconditioner context 664 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 665 PC_MG_FULL, PC_MG_KASKADE 666 667 Options Database Key: 668 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 669 additive, full, kaskade 670 671 Level: advanced 672 673 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 674 675 .seealso: PCMGSetLevels() 676 @*/ 677 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form) 678 { 679 PC_MG *mg = (PC_MG*)pc->data; 680 681 PetscFunctionBegin; 682 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 683 mg->am = form; 684 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 685 else pc->ops->applyrichardson = 0; 686 PetscFunctionReturn(0); 687 } 688 689 #undef __FUNCT__ 690 #define __FUNCT__ "PCMGSetCycleType" 691 /*@ 692 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 693 complicated cycling. 694 695 Collective on PC 696 697 Input Parameters: 698 + pc - the multigrid context 699 - PC_MG_CYCLE_V or PC_MG_CYCLE_W 700 701 Options Database Key: 702 $ -pc_mg_cycle_type v or w 703 704 Level: advanced 705 706 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 707 708 .seealso: PCMGSetCycleTypeOnLevel() 709 @*/ 710 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n) 711 { 712 PC_MG *mg = (PC_MG*)pc->data; 713 PC_MG_Levels **mglevels = mg->levels; 714 PetscInt i,levels; 715 716 PetscFunctionBegin; 717 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 718 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 719 levels = mglevels[0]->levels; 720 721 for (i=0; i<levels; i++) { 722 mglevels[i]->cycles = n; 723 } 724 PetscFunctionReturn(0); 725 } 726 727 #undef __FUNCT__ 728 #define __FUNCT__ "PCMGMultiplicativeSetCycles" 729 /*@ 730 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 731 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 732 733 Collective on PC 734 735 Input Parameters: 736 + pc - the multigrid context 737 - n - number of cycles (default is 1) 738 739 Options Database Key: 740 $ -pc_mg_multiplicative_cycles n 741 742 Level: advanced 743 744 Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 745 746 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 747 748 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 749 @*/ 750 PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 751 { 752 PC_MG *mg = (PC_MG*)pc->data; 753 PC_MG_Levels **mglevels = mg->levels; 754 PetscInt i,levels; 755 756 PetscFunctionBegin; 757 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 758 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 759 levels = mglevels[0]->levels; 760 761 for (i=0; i<levels; i++) { 762 mg->cyclesperpcapply = n; 763 } 764 PetscFunctionReturn(0); 765 } 766 767 #undef __FUNCT__ 768 #define __FUNCT__ "PCMGSetGalerkin" 769 /*@ 770 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 771 finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 772 773 Collective on PC 774 775 Input Parameters: 776 . pc - the multigrid context 777 778 Options Database Key: 779 $ -pc_mg_galerkin 780 781 Level: intermediate 782 783 .keywords: MG, set, Galerkin 784 785 .seealso: PCMGGetGalerkin() 786 787 @*/ 788 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 789 { 790 PC_MG *mg = (PC_MG*)pc->data; 791 PC_MG_Levels **mglevels = mg->levels; 792 PetscInt i,levels; 793 794 PetscFunctionBegin; 795 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 796 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 797 levels = mglevels[0]->levels; 798 799 for (i=0; i<levels; i++) { 800 mglevels[i]->galerkin = PETSC_TRUE; 801 } 802 PetscFunctionReturn(0); 803 } 804 805 #undef __FUNCT__ 806 #define __FUNCT__ "PCMGGetGalerkin" 807 /*@ 808 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 809 A_i-1 = r_i * A_i * r_i^t 810 811 Not Collective 812 813 Input Parameter: 814 . pc - the multigrid context 815 816 Output Parameter: 817 . gelerkin - PETSC_TRUE or PETSC_FALSE 818 819 Options Database Key: 820 $ -pc_mg_galerkin 821 822 Level: intermediate 823 824 .keywords: MG, set, Galerkin 825 826 .seealso: PCMGSetGalerkin() 827 828 @*/ 829 PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 830 { 831 PC_MG *mg = (PC_MG*)pc->data; 832 PC_MG_Levels **mglevels = mg->levels; 833 834 PetscFunctionBegin; 835 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 836 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 837 *galerkin = mglevels[0]->galerkin; 838 PetscFunctionReturn(0); 839 } 840 841 #undef __FUNCT__ 842 #define __FUNCT__ "PCMGSetNumberSmoothDown" 843 /*@ 844 PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 845 use on all levels. Use PCMGGetSmootherDown() to set different 846 pre-smoothing steps on different levels. 847 848 Collective on PC 849 850 Input Parameters: 851 + mg - the multigrid context 852 - n - the number of smoothing steps 853 854 Options Database Key: 855 . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 856 857 Level: advanced 858 859 .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 860 861 .seealso: PCMGSetNumberSmoothUp() 862 @*/ 863 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n) 864 { 865 PC_MG *mg = (PC_MG*)pc->data; 866 PC_MG_Levels **mglevels = mg->levels; 867 PetscErrorCode ierr; 868 PetscInt i,levels; 869 870 PetscFunctionBegin; 871 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 872 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 873 levels = mglevels[0]->levels; 874 875 for (i=1; i<levels; i++) { 876 /* make sure smoother up and down are different */ 877 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 878 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 879 mg->default_smoothd = n; 880 } 881 PetscFunctionReturn(0); 882 } 883 884 #undef __FUNCT__ 885 #define __FUNCT__ "PCMGSetNumberSmoothUp" 886 /*@ 887 PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 888 on all levels. Use PCMGGetSmootherUp() to set different numbers of 889 post-smoothing steps on different levels. 890 891 Collective on PC 892 893 Input Parameters: 894 + mg - the multigrid context 895 - n - the number of smoothing steps 896 897 Options Database Key: 898 . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 899 900 Level: advanced 901 902 Note: this does not set a value on the coarsest grid, since we assume that 903 there is no separate smooth up on the coarsest grid. 904 905 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 906 907 .seealso: PCMGSetNumberSmoothDown() 908 @*/ 909 PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 910 { 911 PC_MG *mg = (PC_MG*)pc->data; 912 PC_MG_Levels **mglevels = mg->levels; 913 PetscErrorCode ierr; 914 PetscInt i,levels; 915 916 PetscFunctionBegin; 917 PetscValidHeaderSpecific(pc,PC_COOKIE,1); 918 if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 919 levels = mglevels[0]->levels; 920 921 for (i=1; i<levels; i++) { 922 /* make sure smoother up and down are different */ 923 ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 924 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 925 mg->default_smoothu = n; 926 } 927 PetscFunctionReturn(0); 928 } 929 930 /* ----------------------------------------------------------------------------------------*/ 931 932 /*MC 933 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 934 information about the coarser grid matrices and restriction/interpolation operators. 935 936 Options Database Keys: 937 + -pc_mg_levels <nlevels> - number of levels including finest 938 . -pc_mg_cycles v or w 939 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 940 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 941 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 942 . -pc_mg_log - log information about time spent on each level of the solver 943 . -pc_mg_monitor - print information on the multigrid convergence 944 . -pc_mg_galerkin - use Galerkin process to compute coarser operators 945 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 946 to the Socket viewer for reading from Matlab. 947 948 Notes: 949 950 Level: intermediate 951 952 Concepts: multigrid/multilevel 953 954 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 955 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 956 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 957 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 958 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 959 M*/ 960 961 EXTERN_C_BEGIN 962 #undef __FUNCT__ 963 #define __FUNCT__ "PCCreate_MG" 964 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 965 { 966 PC_MG *mg; 967 PetscErrorCode ierr; 968 969 PetscFunctionBegin; 970 ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 971 pc->data = (void*)mg; 972 mg->nlevels = -1; 973 974 pc->ops->apply = PCApply_MG; 975 pc->ops->setup = PCSetUp_MG; 976 pc->ops->destroy = PCDestroy_MG; 977 pc->ops->setfromoptions = PCSetFromOptions_MG; 978 pc->ops->view = PCView_MG; 979 PetscFunctionReturn(0); 980 } 981 EXTERN_C_END 982