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