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