1 2 /* 3 Defines the multigrid preconditioner interface. 4 */ 5 #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 6 #include <petscdm.h> 7 PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*); 8 9 PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 10 { 11 PC_MG *mg = (PC_MG*)pc->data; 12 PC_MG_Levels *mgc,*mglevels = *mglevelsin; 13 PetscErrorCode ierr; 14 PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 15 PC subpc; 16 PCFailedReason pcreason; 17 18 PetscFunctionBegin; 19 if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 20 ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 21 ierr = KSPGetPC(mglevels->smoothd,&subpc);CHKERRQ(ierr); 22 ierr = PCGetSetUpFailedReason(subpc,&pcreason);CHKERRQ(ierr); 23 if (pcreason) { 24 pc->failedreason = PC_SUBPC_ERROR; 25 } 26 if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 27 if (mglevels->level) { /* not the coarsest grid */ 28 if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 29 ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 30 if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 31 32 /* if on finest level and have convergence criteria set */ 33 if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 34 PetscReal rnorm; 35 ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 36 if (rnorm <= mg->ttol) { 37 if (rnorm < mg->abstol) { 38 *reason = PCRICHARDSON_CONVERGED_ATOL; 39 ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr); 40 } else { 41 *reason = PCRICHARDSON_CONVERGED_RTOL; 42 ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);CHKERRQ(ierr); 43 } 44 PetscFunctionReturn(0); 45 } 46 } 47 48 mgc = *(mglevelsin - 1); 49 if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 50 ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 51 if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 52 ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 53 while (cycles--) { 54 ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 55 } 56 if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 57 ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 58 if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 59 if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 60 ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 61 if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 62 } 63 PetscFunctionReturn(0); 64 } 65 66 static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscBool zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason) 67 { 68 PC_MG *mg = (PC_MG*)pc->data; 69 PC_MG_Levels **mglevels = mg->levels; 70 PetscErrorCode ierr; 71 PC tpc; 72 PetscBool changeu,changed; 73 PetscInt levels = mglevels[0]->levels,i; 74 75 PetscFunctionBegin; 76 /* When the DM is supplying the matrix then it will not exist until here */ 77 for (i=0; i<levels; i++) { 78 if (!mglevels[i]->A) { 79 ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 80 ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 81 } 82 } 83 84 ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 85 ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 86 ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 87 ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 88 if (!changed && !changeu) { 89 ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 90 mglevels[levels-1]->b = b; 91 } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 92 if (!mglevels[levels-1]->b) { 93 Vec *vec; 94 95 ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 96 mglevels[levels-1]->b = *vec; 97 ierr = PetscFree(vec);CHKERRQ(ierr); 98 } 99 ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 100 } 101 mglevels[levels-1]->x = x; 102 103 mg->rtol = rtol; 104 mg->abstol = abstol; 105 mg->dtol = dtol; 106 if (rtol) { 107 /* compute initial residual norm for relative convergence test */ 108 PetscReal rnorm; 109 if (zeroguess) { 110 ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 111 } else { 112 ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 113 ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 114 } 115 mg->ttol = PetscMax(rtol*rnorm,abstol); 116 } else if (abstol) mg->ttol = abstol; 117 else mg->ttol = 0.0; 118 119 /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 120 stop prematurely due to small residual */ 121 for (i=1; i<levels; i++) { 122 ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 123 if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 124 /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 125 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 126 ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 127 } 128 } 129 130 *reason = (PCRichardsonConvergedReason)0; 131 for (i=0; i<its; i++) { 132 ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 133 if (*reason) break; 134 } 135 if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 136 *outits = i; 137 if (!changed && !changeu) mglevels[levels-1]->b = NULL; 138 PetscFunctionReturn(0); 139 } 140 141 PetscErrorCode PCReset_MG(PC pc) 142 { 143 PC_MG *mg = (PC_MG*)pc->data; 144 PC_MG_Levels **mglevels = mg->levels; 145 PetscErrorCode ierr; 146 PetscInt i,n; 147 148 PetscFunctionBegin; 149 if (mglevels) { 150 n = mglevels[0]->levels; 151 for (i=0; i<n-1; i++) { 152 ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 153 ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 154 ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 155 ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 156 ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 157 ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr); 158 ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 159 } 160 /* this is not null only if the smoother on the finest level 161 changes the rhs during PreSolve */ 162 ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr); 163 164 for (i=0; i<n; i++) { 165 ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 166 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 167 ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 168 } 169 ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 170 } 171 } 172 PetscFunctionReturn(0); 173 } 174 175 PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms) 176 { 177 PetscErrorCode ierr; 178 PC_MG *mg = (PC_MG*)pc->data; 179 MPI_Comm comm; 180 PC_MG_Levels **mglevels = mg->levels; 181 PCMGType mgtype = mg->am; 182 PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V; 183 PetscInt i; 184 PetscMPIInt size; 185 const char *prefix; 186 PC ipc; 187 PetscInt n; 188 189 PetscFunctionBegin; 190 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 191 PetscValidLogicalCollectiveInt(pc,levels,2); 192 if (mg->nlevels == levels) PetscFunctionReturn(0); 193 ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 194 if (mglevels) { 195 mgctype = mglevels[0]->cycles; 196 /* changing the number of levels so free up the previous stuff */ 197 ierr = PCReset_MG(pc);CHKERRQ(ierr); 198 n = mglevels[0]->levels; 199 for (i=0; i<n; i++) { 200 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 201 ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 202 } 203 ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 204 ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 205 } 206 ierr = PetscFree(mg->levels);CHKERRQ(ierr); 207 } 208 209 mg->nlevels = levels; 210 211 ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr); 212 ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 213 214 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 215 216 mg->stageApply = 0; 217 for (i=0; i<levels; i++) { 218 ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr); 219 220 mglevels[i]->level = i; 221 mglevels[i]->levels = levels; 222 mglevels[i]->cycles = mgctype; 223 mg->default_smoothu = 2; 224 mg->default_smoothd = 2; 225 mglevels[i]->eventsmoothsetup = 0; 226 mglevels[i]->eventsmoothsolve = 0; 227 mglevels[i]->eventresidual = 0; 228 mglevels[i]->eventinterprestrict = 0; 229 230 if (comms) comm = comms[i]; 231 ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 232 ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr); 233 ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 234 ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 235 ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 236 if (i || levels == 1) { 237 char tprefix[128]; 238 239 ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 240 ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 241 ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 242 ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 243 ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 244 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 245 246 sprintf(tprefix,"mg_levels_%d_",(int)i); 247 ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 248 } else { 249 ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 250 251 /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 252 ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 253 ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 254 ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 255 if (size > 1) { 256 ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 257 } else { 258 ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 259 } 260 ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 261 } 262 ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr); 263 264 mglevels[i]->smoothu = mglevels[i]->smoothd; 265 mg->rtol = 0.0; 266 mg->abstol = 0.0; 267 mg->dtol = 0.0; 268 mg->ttol = 0.0; 269 mg->cyclesperpcapply = 1; 270 } 271 mg->levels = mglevels; 272 ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 273 PetscFunctionReturn(0); 274 } 275 276 /*@C 277 PCMGSetLevels - Sets the number of levels to use with MG. 278 Must be called before any other MG routine. 279 280 Logically Collective on PC 281 282 Input Parameters: 283 + pc - the preconditioner context 284 . levels - the number of levels 285 - comms - optional communicators for each level; this is to allow solving the coarser problems 286 on smaller sets of processors. 287 288 Level: intermediate 289 290 Notes: 291 If the number of levels is one then the multigrid uses the -mg_levels prefix 292 for setting the level options rather than the -mg_coarse prefix. 293 294 .keywords: MG, set, levels, multigrid 295 296 .seealso: PCMGSetType(), PCMGGetLevels() 297 @*/ 298 PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 299 { 300 PetscErrorCode ierr; 301 302 PetscFunctionBegin; 303 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 304 if (comms) PetscValidPointer(comms,3); 305 ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr); 306 PetscFunctionReturn(0); 307 } 308 309 310 PetscErrorCode PCDestroy_MG(PC pc) 311 { 312 PetscErrorCode ierr; 313 PC_MG *mg = (PC_MG*)pc->data; 314 PC_MG_Levels **mglevels = mg->levels; 315 PetscInt i,n; 316 317 PetscFunctionBegin; 318 ierr = PCReset_MG(pc);CHKERRQ(ierr); 319 if (mglevels) { 320 n = mglevels[0]->levels; 321 for (i=0; i<n; i++) { 322 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 323 ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 324 } 325 ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 326 ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 327 } 328 ierr = PetscFree(mg->levels);CHKERRQ(ierr); 329 } 330 ierr = PetscFree(pc->data);CHKERRQ(ierr); 331 PetscFunctionReturn(0); 332 } 333 334 335 336 extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 337 extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 338 extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 339 340 /* 341 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 342 or full cycle of multigrid. 343 344 Note: 345 A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 346 */ 347 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 348 { 349 PC_MG *mg = (PC_MG*)pc->data; 350 PC_MG_Levels **mglevels = mg->levels; 351 PetscErrorCode ierr; 352 PC tpc; 353 PetscInt levels = mglevels[0]->levels,i; 354 PetscBool changeu,changed; 355 356 PetscFunctionBegin; 357 if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 358 /* When the DM is supplying the matrix then it will not exist until here */ 359 for (i=0; i<levels; i++) { 360 if (!mglevels[i]->A) { 361 ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 362 ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 363 } 364 } 365 366 ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 367 ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 368 ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 369 ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 370 if (!changeu && !changed) { 371 ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 372 mglevels[levels-1]->b = b; 373 } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 374 if (!mglevels[levels-1]->b) { 375 Vec *vec; 376 377 ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 378 mglevels[levels-1]->b = *vec; 379 ierr = PetscFree(vec);CHKERRQ(ierr); 380 } 381 ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 382 } 383 mglevels[levels-1]->x = x; 384 385 if (mg->am == PC_MG_MULTIPLICATIVE) { 386 ierr = VecSet(x,0.0);CHKERRQ(ierr); 387 for (i=0; i<mg->cyclesperpcapply; i++) { 388 ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr); 389 } 390 } else if (mg->am == PC_MG_ADDITIVE) { 391 ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 392 } else if (mg->am == PC_MG_KASKADE) { 393 ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 394 } else { 395 ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 396 } 397 if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 398 if (!changeu && !changed) mglevels[levels-1]->b = NULL; 399 PetscFunctionReturn(0); 400 } 401 402 403 PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 404 { 405 PetscErrorCode ierr; 406 PetscInt levels,cycles; 407 PetscBool flg; 408 PC_MG *mg = (PC_MG*)pc->data; 409 PC_MG_Levels **mglevels; 410 PCMGType mgtype; 411 PCMGCycleType mgctype; 412 PCMGGalerkinType gtype; 413 414 PetscFunctionBegin; 415 levels = PetscMax(mg->nlevels,1); 416 ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr); 417 ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 418 if (!flg && !mg->levels && pc->dm) { 419 ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 420 levels++; 421 mg->usedmfornumberoflevels = PETSC_TRUE; 422 } 423 ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 424 mglevels = mg->levels; 425 426 mgctype = (PCMGCycleType) mglevels[0]->cycles; 427 ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 428 if (flg) { 429 ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 430 } 431 gtype = mg->galerkin; 432 ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg);CHKERRQ(ierr); 433 if (flg) { 434 ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr); 435 } 436 flg = PETSC_FALSE; 437 ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr); 438 if (flg) { 439 ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr); 440 } 441 mgtype = mg->am; 442 ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 443 if (flg) { 444 ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 445 } 446 if (mg->am == PC_MG_MULTIPLICATIVE) { 447 ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 448 if (flg) { 449 ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 450 } 451 } 452 flg = PETSC_FALSE; 453 ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 454 if (flg) { 455 PetscInt i; 456 char eventname[128]; 457 458 levels = mglevels[0]->levels; 459 for (i=0; i<levels; i++) { 460 sprintf(eventname,"MGSetup Level %d",(int)i); 461 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 462 sprintf(eventname,"MGSmooth Level %d",(int)i); 463 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 464 if (i) { 465 sprintf(eventname,"MGResid Level %d",(int)i); 466 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 467 sprintf(eventname,"MGInterp Level %d",(int)i); 468 ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 469 } 470 } 471 472 #if defined(PETSC_USE_LOG) 473 { 474 const char *sname = "MG Apply"; 475 PetscStageLog stageLog; 476 PetscInt st; 477 478 ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 479 for (st = 0; st < stageLog->numStages; ++st) { 480 PetscBool same; 481 482 ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 483 if (same) mg->stageApply = st; 484 } 485 if (!mg->stageApply) { 486 ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 487 } 488 } 489 #endif 490 } 491 ierr = PetscOptionsTail();CHKERRQ(ierr); 492 PetscFunctionReturn(0); 493 } 494 495 const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 496 const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 497 const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0}; 498 499 #include <petscdraw.h> 500 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 501 { 502 PC_MG *mg = (PC_MG*)pc->data; 503 PC_MG_Levels **mglevels = mg->levels; 504 PetscErrorCode ierr; 505 PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 506 PetscBool iascii,isbinary,isdraw; 507 508 PetscFunctionBegin; 509 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 510 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 511 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 512 if (iascii) { 513 const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 514 ierr = PetscViewerASCIIPrintf(viewer," type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 515 if (mg->am == PC_MG_MULTIPLICATIVE) { 516 ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 517 } 518 if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 519 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 520 } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 521 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr); 522 } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 523 ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr); 524 } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 525 ierr = PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr); 526 } else { 527 ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 528 } 529 if (mg->view){ 530 ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr); 531 } 532 for (i=0; i<levels; i++) { 533 if (!i) { 534 ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 535 } else { 536 ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 537 } 538 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 539 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 540 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 541 if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 542 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 543 } else if (i) { 544 ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 545 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 546 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 547 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 548 } 549 } 550 } else if (isbinary) { 551 for (i=levels-1; i>=0; i--) { 552 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 553 if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 554 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 555 } 556 } 557 } else if (isdraw) { 558 PetscDraw draw; 559 PetscReal x,w,y,bottom,th; 560 ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 561 ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 562 ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 563 bottom = y - th; 564 for (i=levels-1; i>=0; i--) { 565 if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 566 ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 567 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 568 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 569 } else { 570 w = 0.5*PetscMin(1.0-x,x); 571 ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 572 ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 573 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 574 ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 575 ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 576 ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 577 } 578 ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 579 bottom -= th; 580 } 581 } 582 PetscFunctionReturn(0); 583 } 584 585 #include <petsc/private/dmimpl.h> 586 #include <petsc/private/kspimpl.h> 587 588 /* 589 Calls setup for the KSP on each level 590 */ 591 PetscErrorCode PCSetUp_MG(PC pc) 592 { 593 PC_MG *mg = (PC_MG*)pc->data; 594 PC_MG_Levels **mglevels = mg->levels; 595 PetscErrorCode ierr; 596 PetscInt i,n; 597 PC cpc; 598 PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 599 Mat dA,dB; 600 Vec tvec; 601 DM *dms; 602 PetscViewer viewer = 0; 603 PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE; 604 605 PetscFunctionBegin; 606 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up"); 607 n = mglevels[0]->levels; 608 /* FIX: Move this to PCSetFromOptions_MG? */ 609 if (mg->usedmfornumberoflevels) { 610 PetscInt levels; 611 ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 612 levels++; 613 if (levels > n) { /* the problem is now being solved on a finer grid */ 614 ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 615 n = levels; 616 ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 617 mglevels = mg->levels; 618 } 619 } 620 ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 621 622 623 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 624 /* so use those from global PC */ 625 /* Is this what we always want? What if user wants to keep old one? */ 626 ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 627 if (opsset) { 628 Mat mmat; 629 ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 630 if (mmat == pc->pmat) opsset = PETSC_FALSE; 631 } 632 633 if (!opsset) { 634 ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 635 if(use_amat){ 636 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); 637 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 638 } 639 else { 640 ierr = PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr); 641 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 642 } 643 } 644 645 for (i=n-1; i>0; i--) { 646 if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 647 missinginterpolate = PETSC_TRUE; 648 continue; 649 } 650 } 651 652 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 653 if (dA == dB) dAeqdB = PETSC_TRUE; 654 if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) { 655 needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 656 } 657 658 659 /* 660 Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 661 Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 662 */ 663 if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 664 /* construct the interpolation from the DMs */ 665 Mat p; 666 Vec rscale; 667 ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 668 dms[n-1] = pc->dm; 669 /* Separately create them so we do not get DMKSP interference between levels */ 670 for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 671 /* 672 Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level. 673 Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options. 674 But it is safe to use -dm_mat_type. 675 676 The mat type should not be hardcoded like this, we need to find a better way. 677 ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr); 678 */ 679 for (i=n-2; i>-1; i--) { 680 DMKSP kdm; 681 PetscBool dmhasrestrict, dmhasinject; 682 ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 683 if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 684 ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 685 /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 686 * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 687 kdm->ops->computerhs = NULL; 688 kdm->rhsctx = NULL; 689 if (!mglevels[i+1]->interpolate) { 690 ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 691 ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 692 if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 693 ierr = VecDestroy(&rscale);CHKERRQ(ierr); 694 ierr = MatDestroy(&p);CHKERRQ(ierr); 695 } 696 ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr); 697 if (dmhasrestrict && !mglevels[i+1]->restrct){ 698 ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr); 699 ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr); 700 ierr = MatDestroy(&p);CHKERRQ(ierr); 701 } 702 ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr); 703 if (dmhasinject && !mglevels[i+1]->inject){ 704 ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr); 705 ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr); 706 ierr = MatDestroy(&p);CHKERRQ(ierr); 707 } 708 } 709 710 for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 711 ierr = PetscFree(dms);CHKERRQ(ierr); 712 } 713 714 if (pc->dm && !pc->setupcalled) { 715 /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 716 ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 717 ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 718 } 719 720 if (mg->galerkin < PC_MG_GALERKIN_NONE) { 721 Mat A,B; 722 PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 723 MatReuse reuse = MAT_INITIAL_MATRIX; 724 725 if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 726 if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE; 727 if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 728 for (i=n-2; i>-1; i--) { 729 if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0"); 730 if (!mglevels[i+1]->interpolate) { 731 ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 732 } 733 if (!mglevels[i+1]->restrct) { 734 ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 735 } 736 if (reuse == MAT_REUSE_MATRIX) { 737 ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr); 738 } 739 if (doA) { 740 ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr); 741 } 742 if (doB) { 743 ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr); 744 } 745 /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 746 if (!doA && dAeqdB) { 747 if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);} 748 A = B; 749 } else if (!doA && reuse == MAT_INITIAL_MATRIX ) { 750 ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr); 751 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 752 } 753 if (!doB && dAeqdB) { 754 if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);} 755 B = A; 756 } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 757 ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 758 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 759 } 760 if (reuse == MAT_INITIAL_MATRIX) { 761 ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr); 762 ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr); 763 ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr); 764 } 765 dA = A; 766 dB = B; 767 } 768 } 769 if (needRestricts && pc->dm && pc->dm->x) { 770 /* need to restrict Jacobian location to coarser meshes for evaluation */ 771 for (i=n-2; i>-1; i--) { 772 Mat R; 773 Vec rscale; 774 if (!mglevels[i]->smoothd->dm->x) { 775 Vec *vecs; 776 ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 777 mglevels[i]->smoothd->dm->x = vecs[0]; 778 ierr = PetscFree(vecs);CHKERRQ(ierr); 779 } 780 ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 781 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 782 ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 783 ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 784 } 785 } 786 if (needRestricts && pc->dm) { 787 for (i=n-2; i>=0; i--) { 788 DM dmfine,dmcoarse; 789 Mat Restrict,Inject; 790 Vec rscale; 791 ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 792 ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 793 ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 794 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 795 ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr); 796 ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 797 } 798 } 799 800 if (!pc->setupcalled) { 801 for (i=0; i<n; i++) { 802 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 803 } 804 for (i=1; i<n; i++) { 805 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 806 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 807 } 808 } 809 /* insure that if either interpolation or restriction is set the other other one is set */ 810 for (i=1; i<n; i++) { 811 ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr); 812 ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr); 813 } 814 for (i=0; i<n-1; i++) { 815 if (!mglevels[i]->b) { 816 Vec *vec; 817 ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 818 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 819 ierr = VecDestroy(vec);CHKERRQ(ierr); 820 ierr = PetscFree(vec);CHKERRQ(ierr); 821 } 822 if (!mglevels[i]->r && i) { 823 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 824 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 825 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 826 } 827 if (!mglevels[i]->x) { 828 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 829 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 830 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 831 } 832 } 833 if (n != 1 && !mglevels[n-1]->r) { 834 /* PCMGSetR() on the finest level if user did not supply it */ 835 Vec *vec; 836 ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 837 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 838 ierr = VecDestroy(vec);CHKERRQ(ierr); 839 ierr = PetscFree(vec);CHKERRQ(ierr); 840 } 841 } 842 843 if (pc->dm) { 844 /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 845 for (i=0; i<n-1; i++) { 846 if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 847 } 848 } 849 850 for (i=1; i<n; i++) { 851 if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 852 /* if doing only down then initial guess is zero */ 853 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 854 } 855 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 856 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 857 if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 858 pc->failedreason = PC_SUBPC_ERROR; 859 } 860 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 861 if (!mglevels[i]->residual) { 862 Mat mat; 863 ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 864 ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 865 } 866 } 867 for (i=1; i<n; i++) { 868 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 869 Mat downmat,downpmat; 870 871 /* check if operators have been set for up, if not use down operators to set them */ 872 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 873 if (!opsset) { 874 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 875 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 876 } 877 878 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 879 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 880 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 881 if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) { 882 pc->failedreason = PC_SUBPC_ERROR; 883 } 884 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 885 } 886 } 887 888 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 889 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 890 if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 891 pc->failedreason = PC_SUBPC_ERROR; 892 } 893 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 894 895 /* 896 Dump the interpolation/restriction matrices plus the 897 Jacobian/stiffness on each level. This allows MATLAB users to 898 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 899 900 Only support one or the other at the same time. 901 */ 902 #if defined(PETSC_USE_SOCKET_VIEWER) 903 ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 904 if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 905 dump = PETSC_FALSE; 906 #endif 907 ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 908 if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 909 910 if (viewer) { 911 for (i=1; i<n; i++) { 912 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 913 } 914 for (i=0; i<n; i++) { 915 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 916 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 917 } 918 } 919 PetscFunctionReturn(0); 920 } 921 922 /* -------------------------------------------------------------------------------------*/ 923 924 PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 925 { 926 PC_MG *mg = (PC_MG *) pc->data; 927 928 PetscFunctionBegin; 929 *levels = mg->nlevels; 930 PetscFunctionReturn(0); 931 } 932 933 /*@ 934 PCMGGetLevels - Gets the number of levels to use with MG. 935 936 Not Collective 937 938 Input Parameter: 939 . pc - the preconditioner context 940 941 Output parameter: 942 . levels - the number of levels 943 944 Level: advanced 945 946 .keywords: MG, get, levels, multigrid 947 948 .seealso: PCMGSetLevels() 949 @*/ 950 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 951 { 952 PetscErrorCode ierr; 953 954 PetscFunctionBegin; 955 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 956 PetscValidIntPointer(levels,2); 957 *levels = 0; 958 ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr); 959 PetscFunctionReturn(0); 960 } 961 962 /*@ 963 PCMGSetType - Determines the form of multigrid to use: 964 multiplicative, additive, full, or the Kaskade algorithm. 965 966 Logically Collective on PC 967 968 Input Parameters: 969 + pc - the preconditioner context 970 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 971 PC_MG_FULL, PC_MG_KASKADE 972 973 Options Database Key: 974 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 975 additive, full, kaskade 976 977 Level: advanced 978 979 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 980 981 .seealso: PCMGSetLevels() 982 @*/ 983 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 984 { 985 PC_MG *mg = (PC_MG*)pc->data; 986 987 PetscFunctionBegin; 988 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 989 PetscValidLogicalCollectiveEnum(pc,form,2); 990 mg->am = form; 991 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 992 else pc->ops->applyrichardson = NULL; 993 PetscFunctionReturn(0); 994 } 995 996 /*@ 997 PCMGGetType - Determines the form of multigrid to use: 998 multiplicative, additive, full, or the Kaskade algorithm. 999 1000 Logically Collective on PC 1001 1002 Input Parameter: 1003 . pc - the preconditioner context 1004 1005 Output Parameter: 1006 . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 1007 1008 1009 Level: advanced 1010 1011 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 1012 1013 .seealso: PCMGSetLevels() 1014 @*/ 1015 PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 1016 { 1017 PC_MG *mg = (PC_MG*)pc->data; 1018 1019 PetscFunctionBegin; 1020 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1021 *type = mg->am; 1022 PetscFunctionReturn(0); 1023 } 1024 1025 /*@ 1026 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 1027 complicated cycling. 1028 1029 Logically Collective on PC 1030 1031 Input Parameters: 1032 + pc - the multigrid context 1033 - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 1034 1035 Options Database Key: 1036 . -pc_mg_cycle_type <v,w> - provide the cycle desired 1037 1038 Level: advanced 1039 1040 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 1041 1042 .seealso: PCMGSetCycleTypeOnLevel() 1043 @*/ 1044 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 1045 { 1046 PC_MG *mg = (PC_MG*)pc->data; 1047 PC_MG_Levels **mglevels = mg->levels; 1048 PetscInt i,levels; 1049 1050 PetscFunctionBegin; 1051 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1052 PetscValidLogicalCollectiveEnum(pc,n,2); 1053 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1054 levels = mglevels[0]->levels; 1055 for (i=0; i<levels; i++) mglevels[i]->cycles = n; 1056 PetscFunctionReturn(0); 1057 } 1058 1059 /*@ 1060 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 1061 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 1062 1063 Logically Collective on PC 1064 1065 Input Parameters: 1066 + pc - the multigrid context 1067 - n - number of cycles (default is 1) 1068 1069 Options Database Key: 1070 . -pc_mg_multiplicative_cycles n 1071 1072 Level: advanced 1073 1074 Notes: 1075 This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 1076 1077 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 1078 1079 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 1080 @*/ 1081 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 1082 { 1083 PC_MG *mg = (PC_MG*)pc->data; 1084 1085 PetscFunctionBegin; 1086 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1087 PetscValidLogicalCollectiveInt(pc,n,2); 1088 mg->cyclesperpcapply = n; 1089 PetscFunctionReturn(0); 1090 } 1091 1092 PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1093 { 1094 PC_MG *mg = (PC_MG*)pc->data; 1095 1096 PetscFunctionBegin; 1097 mg->galerkin = use; 1098 PetscFunctionReturn(0); 1099 } 1100 1101 /*@ 1102 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 1103 finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1104 1105 Logically Collective on PC 1106 1107 Input Parameters: 1108 + pc - the multigrid context 1109 - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1110 1111 Options Database Key: 1112 . -pc_mg_galerkin <both,pmat,mat,none> 1113 1114 Level: intermediate 1115 1116 Notes: 1117 Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 1118 use the PCMG construction of the coarser grids. 1119 1120 .keywords: MG, set, Galerkin 1121 1122 .seealso: PCMGGetGalerkin(), PCMGGalerkinType 1123 1124 @*/ 1125 PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1126 { 1127 PetscErrorCode ierr; 1128 1129 PetscFunctionBegin; 1130 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1131 ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr); 1132 PetscFunctionReturn(0); 1133 } 1134 1135 /*@ 1136 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 1137 A_i-1 = r_i * A_i * p_i 1138 1139 Not Collective 1140 1141 Input Parameter: 1142 . pc - the multigrid context 1143 1144 Output Parameter: 1145 . galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL 1146 1147 Level: intermediate 1148 1149 .keywords: MG, set, Galerkin 1150 1151 .seealso: PCMGSetGalerkin(), PCMGGalerkinType 1152 1153 @*/ 1154 PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 1155 { 1156 PC_MG *mg = (PC_MG*)pc->data; 1157 1158 PetscFunctionBegin; 1159 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1160 *galerkin = mg->galerkin; 1161 PetscFunctionReturn(0); 1162 } 1163 1164 /*@ 1165 PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1166 on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1167 pre- and post-smoothing steps. 1168 1169 Logically Collective on PC 1170 1171 Input Parameters: 1172 + mg - the multigrid context 1173 - n - the number of smoothing steps 1174 1175 Options Database Key: 1176 + -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 1177 1178 Level: advanced 1179 1180 Notes: 1181 this does not set a value on the coarsest grid, since we assume that 1182 there is no separate smooth up on the coarsest grid. 1183 1184 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1185 1186 .seealso: PCMGSetDistinctSmoothUp() 1187 @*/ 1188 PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 1189 { 1190 PC_MG *mg = (PC_MG*)pc->data; 1191 PC_MG_Levels **mglevels = mg->levels; 1192 PetscErrorCode ierr; 1193 PetscInt i,levels; 1194 1195 PetscFunctionBegin; 1196 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1197 PetscValidLogicalCollectiveInt(pc,n,2); 1198 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1199 levels = mglevels[0]->levels; 1200 1201 for (i=1; i<levels; i++) { 1202 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1203 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1204 mg->default_smoothu = n; 1205 mg->default_smoothd = n; 1206 } 1207 PetscFunctionReturn(0); 1208 } 1209 1210 /*@ 1211 PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels 1212 and adds the suffix _up to the options name 1213 1214 Logically Collective on PC 1215 1216 Input Parameters: 1217 . pc - the preconditioner context 1218 1219 Options Database Key: 1220 . -pc_mg_distinct_smoothup 1221 1222 Level: advanced 1223 1224 Notes: 1225 this does not set a value on the coarsest grid, since we assume that 1226 there is no separate smooth up on the coarsest grid. 1227 1228 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1229 1230 .seealso: PCMGSetNumberSmooth() 1231 @*/ 1232 PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1233 { 1234 PC_MG *mg = (PC_MG*)pc->data; 1235 PC_MG_Levels **mglevels = mg->levels; 1236 PetscErrorCode ierr; 1237 PetscInt i,levels; 1238 KSP subksp; 1239 1240 PetscFunctionBegin; 1241 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1242 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1243 levels = mglevels[0]->levels; 1244 1245 for (i=1; i<levels; i++) { 1246 const char *prefix = NULL; 1247 /* make sure smoother up and down are different */ 1248 ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr); 1249 ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr); 1250 ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr); 1251 ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr); 1252 } 1253 PetscFunctionReturn(0); 1254 } 1255 1256 /* ----------------------------------------------------------------------------------------*/ 1257 1258 /*MC 1259 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1260 information about the coarser grid matrices and restriction/interpolation operators. 1261 1262 Options Database Keys: 1263 + -pc_mg_levels <nlevels> - number of levels including finest 1264 . -pc_mg_cycle_type <v,w> - provide the cycle desired 1265 . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 1266 . -pc_mg_log - log information about time spent on each level of the solver 1267 . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 1268 . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 1269 . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 1270 . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1271 to the Socket viewer for reading from MATLAB. 1272 - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 1273 to the binary output file called binaryoutput 1274 1275 Notes: 1276 If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method 1277 1278 When run with a single level the smoother options are used on that level NOT the coarse grid solver options 1279 1280 When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 1281 is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 1282 (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 1283 residual is computed at the end of each cycle. 1284 1285 Level: intermediate 1286 1287 Concepts: multigrid/multilevel 1288 1289 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1290 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), 1291 PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1292 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1293 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1294 M*/ 1295 1296 PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 1297 { 1298 PC_MG *mg; 1299 PetscErrorCode ierr; 1300 1301 PetscFunctionBegin; 1302 ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1303 pc->data = (void*)mg; 1304 mg->nlevels = -1; 1305 mg->am = PC_MG_MULTIPLICATIVE; 1306 mg->galerkin = PC_MG_GALERKIN_NONE; 1307 1308 pc->useAmat = PETSC_TRUE; 1309 1310 pc->ops->apply = PCApply_MG; 1311 pc->ops->setup = PCSetUp_MG; 1312 pc->ops->reset = PCReset_MG; 1313 pc->ops->destroy = PCDestroy_MG; 1314 pc->ops->setfromoptions = PCSetFromOptions_MG; 1315 pc->ops->view = PCView_MG; 1316 1317 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 1318 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr); 1319 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr); 1320 PetscFunctionReturn(0); 1321 } 1322