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