1 2 /* 3 Defines the multigrid preconditioner interface. 4 */ 5 #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 6 #include <petsc/private/kspimpl.h> 7 #include <petscdm.h> 8 PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*); 9 10 /* 11 Contains the list of registered coarse space construction routines 12 */ 13 PetscFunctionList PCMGCoarseList = NULL; 14 15 PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PetscBool transpose,PetscBool matapp,PCRichardsonConvergedReason *reason) 16 { 17 PC_MG *mg = (PC_MG*)pc->data; 18 PC_MG_Levels *mgc,*mglevels = *mglevelsin; 19 PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 20 21 PetscFunctionBegin; 22 if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0)); 23 if (!transpose) { 24 if (matapp) { 25 PetscCall(KSPMatSolve(mglevels->smoothd,mglevels->B,mglevels->X)); /* pre-smooth */ 26 PetscCall(KSPCheckSolve(mglevels->smoothd,pc,NULL)); 27 } else { 28 PetscCall(KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x)); /* pre-smooth */ 29 PetscCall(KSPCheckSolve(mglevels->smoothd,pc,mglevels->x)); 30 } 31 } else { 32 PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 33 PetscCall(KSPSolveTranspose(mglevels->smoothu,mglevels->b,mglevels->x)); /* transpose of post-smooth */ 34 PetscCall(KSPCheckSolve(mglevels->smoothu,pc,mglevels->x)); 35 } 36 if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0)); 37 if (mglevels->level) { /* not the coarsest grid */ 38 if (mglevels->eventresidual) PetscCall(PetscLogEventBegin(mglevels->eventresidual,0,0,0,0)); 39 if (matapp && !mglevels->R) { 40 PetscCall(MatDuplicate(mglevels->B,MAT_DO_NOT_COPY_VALUES,&mglevels->R)); 41 } 42 if (!transpose) { 43 if (matapp) PetscCall((*mglevels->matresidual)(mglevels->A,mglevels->B,mglevels->X,mglevels->R)); 44 else PetscCall((*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r)); 45 } else { 46 if (matapp) PetscCall((*mglevels->matresidualtranspose)(mglevels->A,mglevels->B,mglevels->X,mglevels->R)); 47 else PetscCall((*mglevels->residualtranspose)(mglevels->A,mglevels->b,mglevels->x,mglevels->r)); 48 } 49 if (mglevels->eventresidual) PetscCall(PetscLogEventEnd(mglevels->eventresidual,0,0,0,0)); 50 51 /* if on finest level and have convergence criteria set */ 52 if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 53 PetscReal rnorm; 54 PetscCall(VecNorm(mglevels->r,NORM_2,&rnorm)); 55 if (rnorm <= mg->ttol) { 56 if (rnorm < mg->abstol) { 57 *reason = PCRICHARDSON_CONVERGED_ATOL; 58 PetscCall(PetscInfo(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol)); 59 } else { 60 *reason = PCRICHARDSON_CONVERGED_RTOL; 61 PetscCall(PetscInfo(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol)); 62 } 63 PetscFunctionReturn(0); 64 } 65 } 66 67 mgc = *(mglevelsin - 1); 68 if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0)); 69 if (!transpose) { 70 if (matapp) PetscCall(MatMatRestrict(mglevels->restrct,mglevels->R,&mgc->B)); 71 else PetscCall(MatRestrict(mglevels->restrct,mglevels->r,mgc->b)); 72 } else { 73 if (matapp) PetscCall(MatMatRestrict(mglevels->interpolate,mglevels->R,&mgc->B)); 74 else PetscCall(MatRestrict(mglevels->interpolate,mglevels->r,mgc->b)); 75 } 76 if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0)); 77 if (matapp) { 78 if (!mgc->X) { 79 PetscCall(MatDuplicate(mgc->B,MAT_DO_NOT_COPY_VALUES,&mgc->X)); 80 } else { 81 PetscCall(MatZeroEntries(mgc->X)); 82 } 83 } else { 84 PetscCall(VecZeroEntries(mgc->x)); 85 } 86 while (cycles--) { 87 PetscCall(PCMGMCycle_Private(pc,mglevelsin-1,transpose,matapp,reason)); 88 } 89 if (mglevels->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0)); 90 if (!transpose) { 91 if (matapp) PetscCall(MatMatInterpolateAdd(mglevels->interpolate,mgc->X,mglevels->X,&mglevels->X)); 92 else PetscCall(MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x)); 93 } else { 94 PetscCall(MatInterpolateAdd(mglevels->restrct,mgc->x,mglevels->x,mglevels->x)); 95 } 96 if (mglevels->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0)); 97 if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0)); 98 if (!transpose) { 99 if (matapp) { 100 PetscCall(KSPMatSolve(mglevels->smoothu,mglevels->B,mglevels->X)); /* post smooth */ 101 PetscCall(KSPCheckSolve(mglevels->smoothu,pc,NULL)); 102 } else { 103 PetscCall(KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x)); /* post smooth */ 104 PetscCall(KSPCheckSolve(mglevels->smoothu,pc,mglevels->x)); 105 } 106 } else { 107 PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 108 PetscCall(KSPSolveTranspose(mglevels->smoothd,mglevels->b,mglevels->x)); /* post smooth */ 109 PetscCall(KSPCheckSolve(mglevels->smoothd,pc,mglevels->x)); 110 } 111 if (mglevels->cr) { 112 PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 113 /* TODO Turn on copy and turn off noisy if we have an exact solution 114 PetscCall(VecCopy(mglevels->x, mglevels->crx)); 115 PetscCall(VecCopy(mglevels->b, mglevels->crb)); */ 116 PetscCall(KSPSetNoisy_Private(mglevels->crx)); 117 PetscCall(KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx)); /* compatible relaxation */ 118 PetscCall(KSPCheckSolve(mglevels->cr,pc,mglevels->crx)); 119 } 120 if (mglevels->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0)); 121 } 122 PetscFunctionReturn(0); 123 } 124 125 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) 126 { 127 PC_MG *mg = (PC_MG*)pc->data; 128 PC_MG_Levels **mglevels = mg->levels; 129 PC tpc; 130 PetscBool changeu,changed; 131 PetscInt levels = mglevels[0]->levels,i; 132 133 PetscFunctionBegin; 134 /* When the DM is supplying the matrix then it will not exist until here */ 135 for (i=0; i<levels; i++) { 136 if (!mglevels[i]->A) { 137 PetscCall(KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL)); 138 PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 139 } 140 } 141 142 PetscCall(KSPGetPC(mglevels[levels-1]->smoothd,&tpc)); 143 PetscCall(PCPreSolveChangeRHS(tpc,&changed)); 144 PetscCall(KSPGetPC(mglevels[levels-1]->smoothu,&tpc)); 145 PetscCall(PCPreSolveChangeRHS(tpc,&changeu)); 146 if (!changed && !changeu) { 147 PetscCall(VecDestroy(&mglevels[levels-1]->b)); 148 mglevels[levels-1]->b = b; 149 } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 150 if (!mglevels[levels-1]->b) { 151 Vec *vec; 152 153 PetscCall(KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL)); 154 mglevels[levels-1]->b = *vec; 155 PetscCall(PetscFree(vec)); 156 } 157 PetscCall(VecCopy(b,mglevels[levels-1]->b)); 158 } 159 mglevels[levels-1]->x = x; 160 161 mg->rtol = rtol; 162 mg->abstol = abstol; 163 mg->dtol = dtol; 164 if (rtol) { 165 /* compute initial residual norm for relative convergence test */ 166 PetscReal rnorm; 167 if (zeroguess) { 168 PetscCall(VecNorm(b,NORM_2,&rnorm)); 169 } else { 170 PetscCall((*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w)); 171 PetscCall(VecNorm(w,NORM_2,&rnorm)); 172 } 173 mg->ttol = PetscMax(rtol*rnorm,abstol); 174 } else if (abstol) mg->ttol = abstol; 175 else mg->ttol = 0.0; 176 177 /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 178 stop prematurely due to small residual */ 179 for (i=1; i<levels; i++) { 180 PetscCall(KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 181 if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 182 /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 183 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE)); 184 PetscCall(KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT)); 185 } 186 } 187 188 *reason = (PCRichardsonConvergedReason)0; 189 for (i=0; i<its; i++) { 190 PetscCall(PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_FALSE,PETSC_FALSE,reason)); 191 if (*reason) break; 192 } 193 if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 194 *outits = i; 195 if (!changed && !changeu) mglevels[levels-1]->b = NULL; 196 PetscFunctionReturn(0); 197 } 198 199 PetscErrorCode PCReset_MG(PC pc) 200 { 201 PC_MG *mg = (PC_MG*)pc->data; 202 PC_MG_Levels **mglevels = mg->levels; 203 PetscInt i,n; 204 205 PetscFunctionBegin; 206 if (mglevels) { 207 n = mglevels[0]->levels; 208 for (i=0; i<n-1; i++) { 209 PetscCall(VecDestroy(&mglevels[i+1]->r)); 210 PetscCall(VecDestroy(&mglevels[i]->b)); 211 PetscCall(VecDestroy(&mglevels[i]->x)); 212 PetscCall(MatDestroy(&mglevels[i+1]->R)); 213 PetscCall(MatDestroy(&mglevels[i]->B)); 214 PetscCall(MatDestroy(&mglevels[i]->X)); 215 PetscCall(VecDestroy(&mglevels[i]->crx)); 216 PetscCall(VecDestroy(&mglevels[i]->crb)); 217 PetscCall(MatDestroy(&mglevels[i+1]->restrct)); 218 PetscCall(MatDestroy(&mglevels[i+1]->interpolate)); 219 PetscCall(MatDestroy(&mglevels[i+1]->inject)); 220 PetscCall(VecDestroy(&mglevels[i+1]->rscale)); 221 } 222 PetscCall(VecDestroy(&mglevels[n-1]->crx)); 223 PetscCall(VecDestroy(&mglevels[n-1]->crb)); 224 /* this is not null only if the smoother on the finest level 225 changes the rhs during PreSolve */ 226 PetscCall(VecDestroy(&mglevels[n-1]->b)); 227 PetscCall(MatDestroy(&mglevels[n-1]->B)); 228 229 for (i=0; i<n; i++) { 230 PetscCall(MatDestroy(&mglevels[i]->coarseSpace)); 231 PetscCall(MatDestroy(&mglevels[i]->A)); 232 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 233 PetscCall(KSPReset(mglevels[i]->smoothd)); 234 } 235 PetscCall(KSPReset(mglevels[i]->smoothu)); 236 if (mglevels[i]->cr) PetscCall(KSPReset(mglevels[i]->cr)); 237 } 238 mg->Nc = 0; 239 } 240 PetscFunctionReturn(0); 241 } 242 243 /* Implementing CR 244 245 We only want to make corrections that ``do not change'' the coarse solution. What we mean by not changing is that if I prolong my coarse solution to the fine grid and then inject that fine solution back to the coarse grid, I get the same answer. Injection is what Brannick calls R. We want the complementary projector to Inj, which we will call S, after Brannick, so that Inj S = 0. Now the orthogonal projector onto the range of Inj^T is 246 247 Inj^T (Inj Inj^T)^{-1} Inj 248 249 and if Inj is a VecScatter, as it is now in PETSc, we have 250 251 Inj^T Inj 252 253 and 254 255 S = I - Inj^T Inj 256 257 since 258 259 Inj S = Inj - (Inj Inj^T) Inj = 0. 260 261 Brannick suggests 262 263 A \to S^T A S \qquad\mathrm{and}\qquad M \to S^T M S 264 265 but I do not think his :math:`S^T S = I` is correct. Our S is an orthogonal projector, so :math:`S^T S = S^2 = S`. We will use 266 267 M^{-1} A \to S M^{-1} A S 268 269 In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left. 270 271 Check: || Inj P - I ||_F < tol 272 Check: In general, Inj Inj^T = I 273 */ 274 275 typedef struct { 276 PC mg; /* The PCMG object */ 277 PetscInt l; /* The multigrid level for this solver */ 278 Mat Inj; /* The injection matrix */ 279 Mat S; /* I - Inj^T Inj */ 280 } CRContext; 281 282 static PetscErrorCode CRSetup_Private(PC pc) 283 { 284 CRContext *ctx; 285 Mat It; 286 287 PetscFunctionBeginUser; 288 PetscCall(PCShellGetContext(pc, &ctx)); 289 PetscCall(PCMGGetInjection(ctx->mg, ctx->l, &It)); 290 PetscCheck(It,PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG"); 291 PetscCall(MatCreateTranspose(It, &ctx->Inj)); 292 PetscCall(MatCreateNormal(ctx->Inj, &ctx->S)); 293 PetscCall(MatScale(ctx->S, -1.0)); 294 PetscCall(MatShift(ctx->S, 1.0)); 295 PetscFunctionReturn(0); 296 } 297 298 static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y) 299 { 300 CRContext *ctx; 301 302 PetscFunctionBeginUser; 303 PetscCall(PCShellGetContext(pc, &ctx)); 304 PetscCall(MatMult(ctx->S, x, y)); 305 PetscFunctionReturn(0); 306 } 307 308 static PetscErrorCode CRDestroy_Private(PC pc) 309 { 310 CRContext *ctx; 311 312 PetscFunctionBeginUser; 313 PetscCall(PCShellGetContext(pc, &ctx)); 314 PetscCall(MatDestroy(&ctx->Inj)); 315 PetscCall(MatDestroy(&ctx->S)); 316 PetscCall(PetscFree(ctx)); 317 PetscCall(PCShellSetContext(pc, NULL)); 318 PetscFunctionReturn(0); 319 } 320 321 static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr) 322 { 323 CRContext *ctx; 324 325 PetscFunctionBeginUser; 326 PetscCall(PCCreate(PetscObjectComm((PetscObject) pc), cr)); 327 PetscCall(PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)")); 328 PetscCall(PetscCalloc1(1, &ctx)); 329 ctx->mg = pc; 330 ctx->l = l; 331 PetscCall(PCSetType(*cr, PCSHELL)); 332 PetscCall(PCShellSetContext(*cr, ctx)); 333 PetscCall(PCShellSetApply(*cr, CRApply_Private)); 334 PetscCall(PCShellSetSetUp(*cr, CRSetup_Private)); 335 PetscCall(PCShellSetDestroy(*cr, CRDestroy_Private)); 336 PetscFunctionReturn(0); 337 } 338 339 PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms) 340 { 341 PC_MG *mg = (PC_MG*)pc->data; 342 MPI_Comm comm; 343 PC_MG_Levels **mglevels = mg->levels; 344 PCMGType mgtype = mg->am; 345 PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V; 346 PetscInt i; 347 PetscMPIInt size; 348 const char *prefix; 349 PC ipc; 350 PetscInt n; 351 352 PetscFunctionBegin; 353 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 354 PetscValidLogicalCollectiveInt(pc,levels,2); 355 if (mg->nlevels == levels) PetscFunctionReturn(0); 356 PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 357 if (mglevels) { 358 mgctype = mglevels[0]->cycles; 359 /* changing the number of levels so free up the previous stuff */ 360 PetscCall(PCReset_MG(pc)); 361 n = mglevels[0]->levels; 362 for (i=0; i<n; i++) { 363 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 364 PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 365 } 366 PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 367 PetscCall(KSPDestroy(&mglevels[i]->cr)); 368 PetscCall(PetscFree(mglevels[i])); 369 } 370 PetscCall(PetscFree(mg->levels)); 371 } 372 373 mg->nlevels = levels; 374 375 PetscCall(PetscMalloc1(levels,&mglevels)); 376 PetscCall(PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)))); 377 378 PetscCall(PCGetOptionsPrefix(pc,&prefix)); 379 380 mg->stageApply = 0; 381 for (i=0; i<levels; i++) { 382 PetscCall(PetscNewLog(pc,&mglevels[i])); 383 384 mglevels[i]->level = i; 385 mglevels[i]->levels = levels; 386 mglevels[i]->cycles = mgctype; 387 mg->default_smoothu = 2; 388 mg->default_smoothd = 2; 389 mglevels[i]->eventsmoothsetup = 0; 390 mglevels[i]->eventsmoothsolve = 0; 391 mglevels[i]->eventresidual = 0; 392 mglevels[i]->eventinterprestrict = 0; 393 394 if (comms) comm = comms[i]; 395 if (comm != MPI_COMM_NULL) { 396 PetscCall(KSPCreate(comm,&mglevels[i]->smoothd)); 397 PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure)); 398 PetscCall(PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i)); 399 PetscCall(KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix)); 400 PetscCall(PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level)); 401 if (i || levels == 1) { 402 char tprefix[128]; 403 404 PetscCall(KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV)); 405 PetscCall(KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL)); 406 PetscCall(KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE)); 407 PetscCall(KSPGetPC(mglevels[i]->smoothd,&ipc)); 408 PetscCall(PCSetType(ipc,PCSOR)); 409 PetscCall(KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd)); 410 411 PetscCall(PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i)); 412 PetscCall(KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix)); 413 } else { 414 PetscCall(KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_")); 415 416 /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 417 PetscCall(KSPSetType(mglevels[0]->smoothd,KSPPREONLY)); 418 PetscCall(KSPGetPC(mglevels[0]->smoothd,&ipc)); 419 PetscCallMPI(MPI_Comm_size(comm,&size)); 420 if (size > 1) { 421 PetscCall(PCSetType(ipc,PCREDUNDANT)); 422 } else { 423 PetscCall(PCSetType(ipc,PCLU)); 424 } 425 PetscCall(PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS)); 426 } 427 PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd)); 428 } 429 mglevels[i]->smoothu = mglevels[i]->smoothd; 430 mg->rtol = 0.0; 431 mg->abstol = 0.0; 432 mg->dtol = 0.0; 433 mg->ttol = 0.0; 434 mg->cyclesperpcapply = 1; 435 } 436 mg->levels = mglevels; 437 PetscCall(PCMGSetType(pc,mgtype)); 438 PetscFunctionReturn(0); 439 } 440 441 /*@C 442 PCMGSetLevels - Sets the number of levels to use with MG. 443 Must be called before any other MG routine. 444 445 Logically Collective on PC 446 447 Input Parameters: 448 + pc - the preconditioner context 449 . levels - the number of levels 450 - comms - optional communicators for each level; this is to allow solving the coarser problems 451 on smaller sets of processes. For processes that are not included in the computation 452 you must pass MPI_COMM_NULL. Use comms = NULL to specify that all processes 453 should participate in each level of problem. 454 455 Level: intermediate 456 457 Notes: 458 If the number of levels is one then the multigrid uses the -mg_levels prefix 459 for setting the level options rather than the -mg_coarse prefix. 460 461 You can free the information in comms after this routine is called. 462 463 The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level 464 are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on 465 the two levels, rank 0 in the original communicator will pass in an array of 2 communicators 466 of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators 467 the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate 468 in the coarse grid solve. 469 470 Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one 471 must take special care in providing the restriction and interpolation operation. We recommend 472 providing these as two step operations; first perform a standard restriction or interpolation on 473 the full number of ranks for that level and then use an MPI call to copy the resulting vector 474 array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both 475 cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and 476 recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries. 477 478 Fortran Notes: 479 Use comms = PETSC_NULL_MPI_COMM as the equivalent of NULL in the C interface. Note PETSC_NULL_MPI_COMM 480 is not MPI_COMM_NULL. It is more like PETSC_NULL_INTEGER, PETSC_NULL_REAL etc. 481 482 .seealso: `PCMGSetType()`, `PCMGGetLevels()` 483 @*/ 484 PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 485 { 486 PetscFunctionBegin; 487 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 488 if (comms) PetscValidPointer(comms,3); 489 PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms)); 490 PetscFunctionReturn(0); 491 } 492 493 PetscErrorCode PCDestroy_MG(PC pc) 494 { 495 PC_MG *mg = (PC_MG*)pc->data; 496 PC_MG_Levels **mglevels = mg->levels; 497 PetscInt i,n; 498 499 PetscFunctionBegin; 500 PetscCall(PCReset_MG(pc)); 501 if (mglevels) { 502 n = mglevels[0]->levels; 503 for (i=0; i<n; i++) { 504 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 505 PetscCall(KSPDestroy(&mglevels[i]->smoothd)); 506 } 507 PetscCall(KSPDestroy(&mglevels[i]->smoothu)); 508 PetscCall(KSPDestroy(&mglevels[i]->cr)); 509 PetscCall(PetscFree(mglevels[i])); 510 } 511 PetscCall(PetscFree(mg->levels)); 512 } 513 PetscCall(PetscFree(pc->data)); 514 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL)); 515 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL)); 516 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",NULL)); 517 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",NULL)); 518 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",NULL)); 519 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL)); 520 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL)); 521 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",NULL)); 522 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",NULL)); 523 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",NULL)); 524 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",NULL)); 525 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCoarseSpaceType_C",NULL)); 526 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCoarseSpaceType_C",NULL)); 527 PetscFunctionReturn(0); 528 } 529 530 /* 531 PCApply_MG - Runs either an additive, multiplicative, Kaskadic 532 or full cycle of multigrid. 533 534 Note: 535 A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 536 */ 537 static PetscErrorCode PCApply_MG_Internal(PC pc,Vec b,Vec x,Mat B,Mat X,PetscBool transpose) 538 { 539 PC_MG *mg = (PC_MG*)pc->data; 540 PC_MG_Levels **mglevels = mg->levels; 541 PC tpc; 542 PetscInt levels = mglevels[0]->levels,i; 543 PetscBool changeu,changed,matapp; 544 545 PetscFunctionBegin; 546 matapp = (PetscBool)(B && X); 547 if (mg->stageApply) PetscCall(PetscLogStagePush(mg->stageApply)); 548 /* When the DM is supplying the matrix then it will not exist until here */ 549 for (i=0; i<levels; i++) { 550 if (!mglevels[i]->A) { 551 PetscCall(KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL)); 552 PetscCall(PetscObjectReference((PetscObject)mglevels[i]->A)); 553 } 554 } 555 556 PetscCall(KSPGetPC(mglevels[levels-1]->smoothd,&tpc)); 557 PetscCall(PCPreSolveChangeRHS(tpc,&changed)); 558 PetscCall(KSPGetPC(mglevels[levels-1]->smoothu,&tpc)); 559 PetscCall(PCPreSolveChangeRHS(tpc,&changeu)); 560 if (!changeu && !changed) { 561 if (matapp) { 562 PetscCall(MatDestroy(&mglevels[levels-1]->B)); 563 mglevels[levels-1]->B = B; 564 } else { 565 PetscCall(VecDestroy(&mglevels[levels-1]->b)); 566 mglevels[levels-1]->b = b; 567 } 568 } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 569 if (matapp) { 570 if (mglevels[levels-1]->B) { 571 PetscInt N1,N2; 572 PetscBool flg; 573 574 PetscCall(MatGetSize(mglevels[levels-1]->B,NULL,&N1)); 575 PetscCall(MatGetSize(B,NULL,&N2)); 576 PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels-1]->B,((PetscObject)B)->type_name,&flg)); 577 if (N1 != N2 || !flg) { 578 PetscCall(MatDestroy(&mglevels[levels-1]->B)); 579 } 580 } 581 if (!mglevels[levels-1]->B) { 582 PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&mglevels[levels-1]->B)); 583 } else { 584 PetscCall(MatCopy(B,mglevels[levels-1]->B,SAME_NONZERO_PATTERN)); 585 } 586 } else { 587 if (!mglevels[levels-1]->b) { 588 Vec *vec; 589 590 PetscCall(KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL)); 591 mglevels[levels-1]->b = *vec; 592 PetscCall(PetscFree(vec)); 593 } 594 PetscCall(VecCopy(b,mglevels[levels-1]->b)); 595 } 596 } 597 if (matapp) { mglevels[levels-1]->X = X; } 598 else { mglevels[levels-1]->x = x; } 599 600 /* If coarser Xs are present, it means we have already block applied the PC at least once 601 Reset operators if sizes/type do no match */ 602 if (matapp && levels > 1 && mglevels[levels-2]->X) { 603 PetscInt Xc,Bc; 604 PetscBool flg; 605 606 PetscCall(MatGetSize(mglevels[levels-2]->X,NULL,&Xc)); 607 PetscCall(MatGetSize(mglevels[levels-1]->B,NULL,&Bc)); 608 PetscCall(PetscObjectTypeCompare((PetscObject)mglevels[levels-2]->X,((PetscObject)mglevels[levels-1]->X)->type_name,&flg)); 609 if (Xc != Bc || !flg) { 610 PetscCall(MatDestroy(&mglevels[levels-1]->R)); 611 for (i=0;i<levels-1;i++) { 612 PetscCall(MatDestroy(&mglevels[i]->R)); 613 PetscCall(MatDestroy(&mglevels[i]->B)); 614 PetscCall(MatDestroy(&mglevels[i]->X)); 615 } 616 } 617 } 618 619 if (mg->am == PC_MG_MULTIPLICATIVE) { 620 if (matapp) PetscCall(MatZeroEntries(X)); 621 else PetscCall(VecZeroEntries(x)); 622 for (i=0; i<mg->cyclesperpcapply; i++) { 623 PetscCall(PCMGMCycle_Private(pc,mglevels+levels-1,transpose,matapp,NULL)); 624 } 625 } else if (mg->am == PC_MG_ADDITIVE) { 626 PetscCall(PCMGACycle_Private(pc,mglevels,transpose,matapp)); 627 } else if (mg->am == PC_MG_KASKADE) { 628 PetscCall(PCMGKCycle_Private(pc,mglevels,transpose,matapp)); 629 } else { 630 PetscCall(PCMGFCycle_Private(pc,mglevels,transpose,matapp)); 631 } 632 if (mg->stageApply) PetscCall(PetscLogStagePop()); 633 if (!changeu && !changed) { 634 if (matapp) { mglevels[levels-1]->B = NULL; } 635 else { mglevels[levels-1]->b = NULL; } 636 } 637 PetscFunctionReturn(0); 638 } 639 640 static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 641 { 642 PetscFunctionBegin; 643 PetscCall(PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_FALSE)); 644 PetscFunctionReturn(0); 645 } 646 647 static PetscErrorCode PCApplyTranspose_MG(PC pc,Vec b,Vec x) 648 { 649 PetscFunctionBegin; 650 PetscCall(PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_TRUE)); 651 PetscFunctionReturn(0); 652 } 653 654 static PetscErrorCode PCMatApply_MG(PC pc,Mat b,Mat x) 655 { 656 PetscFunctionBegin; 657 PetscCall(PCApply_MG_Internal(pc,NULL,NULL,b,x,PETSC_FALSE)); 658 PetscFunctionReturn(0); 659 } 660 661 PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 662 { 663 PetscInt levels,cycles; 664 PetscBool flg, flg2; 665 PC_MG *mg = (PC_MG*)pc->data; 666 PC_MG_Levels **mglevels; 667 PCMGType mgtype; 668 PCMGCycleType mgctype; 669 PCMGGalerkinType gtype; 670 PCMGCoarseSpaceType coarseSpaceType; 671 672 PetscFunctionBegin; 673 levels = PetscMax(mg->nlevels,1); 674 PetscOptionsHeadBegin(PetscOptionsObject,"Multigrid options"); 675 PetscCall(PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg)); 676 if (!flg && !mg->levels && pc->dm) { 677 PetscCall(DMGetRefineLevel(pc->dm,&levels)); 678 levels++; 679 mg->usedmfornumberoflevels = PETSC_TRUE; 680 } 681 PetscCall(PCMGSetLevels(pc,levels,NULL)); 682 mglevels = mg->levels; 683 684 mgctype = (PCMGCycleType) mglevels[0]->cycles; 685 PetscCall(PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg)); 686 if (flg) { 687 PetscCall(PCMGSetCycleType(pc,mgctype)); 688 } 689 gtype = mg->galerkin; 690 PetscCall(PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg)); 691 if (flg) { 692 PetscCall(PCMGSetGalerkin(pc,gtype)); 693 } 694 coarseSpaceType = mg->coarseSpaceType; 695 PetscCall(PetscOptionsEnum("-pc_mg_adapt_interp_coarse_space","Type of adaptive coarse space: none, polynomial, harmonic, eigenvector, generalized_eigenvector, gdsw","PCMGSetAdaptCoarseSpaceType",PCMGCoarseSpaceTypes,(PetscEnum)coarseSpaceType,(PetscEnum*)&coarseSpaceType,&flg)); 696 if (flg) PetscCall(PCMGSetAdaptCoarseSpaceType(pc,coarseSpaceType)); 697 PetscCall(PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg)); 698 PetscCall(PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg)); 699 flg2 = PETSC_FALSE; 700 PetscCall(PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg)); 701 if (flg) PetscCall(PCMGSetAdaptCR(pc, flg2)); 702 flg = PETSC_FALSE; 703 PetscCall(PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL)); 704 if (flg) { 705 PetscCall(PCMGSetDistinctSmoothUp(pc)); 706 } 707 mgtype = mg->am; 708 PetscCall(PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg)); 709 if (flg) { 710 PetscCall(PCMGSetType(pc,mgtype)); 711 } 712 if (mg->am == PC_MG_MULTIPLICATIVE) { 713 PetscCall(PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg)); 714 if (flg) { 715 PetscCall(PCMGMultiplicativeSetCycles(pc,cycles)); 716 } 717 } 718 flg = PETSC_FALSE; 719 PetscCall(PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL)); 720 if (flg) { 721 PetscInt i; 722 char eventname[128]; 723 724 levels = mglevels[0]->levels; 725 for (i=0; i<levels; i++) { 726 sprintf(eventname,"MGSetup Level %d",(int)i); 727 PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup)); 728 sprintf(eventname,"MGSmooth Level %d",(int)i); 729 PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve)); 730 if (i) { 731 sprintf(eventname,"MGResid Level %d",(int)i); 732 PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual)); 733 sprintf(eventname,"MGInterp Level %d",(int)i); 734 PetscCall(PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict)); 735 } 736 } 737 738 #if defined(PETSC_USE_LOG) 739 { 740 const char *sname = "MG Apply"; 741 PetscStageLog stageLog; 742 PetscInt st; 743 744 PetscCall(PetscLogGetStageLog(&stageLog)); 745 for (st = 0; st < stageLog->numStages; ++st) { 746 PetscBool same; 747 748 PetscCall(PetscStrcmp(stageLog->stageInfo[st].name, sname, &same)); 749 if (same) mg->stageApply = st; 750 } 751 if (!mg->stageApply) { 752 PetscCall(PetscLogStageRegister(sname, &mg->stageApply)); 753 } 754 } 755 #endif 756 } 757 PetscOptionsHeadEnd(); 758 /* Check option consistency */ 759 PetscCall(PCMGGetGalerkin(pc, >ype)); 760 PetscCall(PCMGGetAdaptInterpolation(pc, &flg)); 761 PetscCheck(!flg || !(gtype >= PC_MG_GALERKIN_NONE),PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator"); 762 PetscFunctionReturn(0); 763 } 764 765 const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL}; 766 const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL}; 767 const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL}; 768 const char *const PCMGCoarseSpaceTypes[] = {"none","polynomial","harmonic","eigenvector","generalized_eigenvector","gdsw","PCMGCoarseSpaceType","PCMG_ADAPT_NONE",NULL}; 769 770 #include <petscdraw.h> 771 PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 772 { 773 PC_MG *mg = (PC_MG*)pc->data; 774 PC_MG_Levels **mglevels = mg->levels; 775 PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 776 PetscBool iascii,isbinary,isdraw; 777 778 PetscFunctionBegin; 779 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 780 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary)); 781 PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 782 if (iascii) { 783 const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 784 PetscCall(PetscViewerASCIIPrintf(viewer," type is %s, levels=%" PetscInt_FMT " cycles=%s\n", PCMGTypes[mg->am],levels,cyclename)); 785 if (mg->am == PC_MG_MULTIPLICATIVE) { 786 PetscCall(PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%" PetscInt_FMT "\n",mg->cyclesperpcapply)); 787 } 788 if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 789 PetscCall(PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n")); 790 } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 791 PetscCall(PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n")); 792 } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 793 PetscCall(PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n")); 794 } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 795 PetscCall(PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n")); 796 } else { 797 PetscCall(PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n")); 798 } 799 if (mg->view) { 800 PetscCall((*mg->view)(pc,viewer)); 801 } 802 for (i=0; i<levels; i++) { 803 if (i) { 804 PetscCall(PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n",i)); 805 } else { 806 PetscCall(PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n",i)); 807 } 808 PetscCall(PetscViewerASCIIPushTab(viewer)); 809 PetscCall(KSPView(mglevels[i]->smoothd,viewer)); 810 PetscCall(PetscViewerASCIIPopTab(viewer)); 811 if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 812 PetscCall(PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n")); 813 } else if (i) { 814 PetscCall(PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n",i)); 815 PetscCall(PetscViewerASCIIPushTab(viewer)); 816 PetscCall(KSPView(mglevels[i]->smoothu,viewer)); 817 PetscCall(PetscViewerASCIIPopTab(viewer)); 818 } 819 if (i && mglevels[i]->cr) { 820 PetscCall(PetscViewerASCIIPrintf(viewer,"CR solver on level %" PetscInt_FMT " -------------------------------\n",i)); 821 PetscCall(PetscViewerASCIIPushTab(viewer)); 822 PetscCall(KSPView(mglevels[i]->cr,viewer)); 823 PetscCall(PetscViewerASCIIPopTab(viewer)); 824 } 825 } 826 } else if (isbinary) { 827 for (i=levels-1; i>=0; i--) { 828 PetscCall(KSPView(mglevels[i]->smoothd,viewer)); 829 if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 830 PetscCall(KSPView(mglevels[i]->smoothu,viewer)); 831 } 832 } 833 } else if (isdraw) { 834 PetscDraw draw; 835 PetscReal x,w,y,bottom,th; 836 PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 837 PetscCall(PetscDrawGetCurrentPoint(draw,&x,&y)); 838 PetscCall(PetscDrawStringGetSize(draw,NULL,&th)); 839 bottom = y - th; 840 for (i=levels-1; i>=0; i--) { 841 if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 842 PetscCall(PetscDrawPushCurrentPoint(draw,x,bottom)); 843 PetscCall(KSPView(mglevels[i]->smoothd,viewer)); 844 PetscCall(PetscDrawPopCurrentPoint(draw)); 845 } else { 846 w = 0.5*PetscMin(1.0-x,x); 847 PetscCall(PetscDrawPushCurrentPoint(draw,x+w,bottom)); 848 PetscCall(KSPView(mglevels[i]->smoothd,viewer)); 849 PetscCall(PetscDrawPopCurrentPoint(draw)); 850 PetscCall(PetscDrawPushCurrentPoint(draw,x-w,bottom)); 851 PetscCall(KSPView(mglevels[i]->smoothu,viewer)); 852 PetscCall(PetscDrawPopCurrentPoint(draw)); 853 } 854 PetscCall(PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL)); 855 bottom -= th; 856 } 857 } 858 PetscFunctionReturn(0); 859 } 860 861 #include <petsc/private/kspimpl.h> 862 863 /* 864 Calls setup for the KSP on each level 865 */ 866 PetscErrorCode PCSetUp_MG(PC pc) 867 { 868 PC_MG *mg = (PC_MG*)pc->data; 869 PC_MG_Levels **mglevels = mg->levels; 870 PetscInt i,n; 871 PC cpc; 872 PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 873 Mat dA,dB; 874 Vec tvec; 875 DM *dms; 876 PetscViewer viewer = NULL; 877 PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE; 878 PetscBool adaptInterpolation = mg->adaptInterpolation; 879 880 PetscFunctionBegin; 881 PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up"); 882 n = mglevels[0]->levels; 883 /* FIX: Move this to PCSetFromOptions_MG? */ 884 if (mg->usedmfornumberoflevels) { 885 PetscInt levels; 886 PetscCall(DMGetRefineLevel(pc->dm,&levels)); 887 levels++; 888 if (levels > n) { /* the problem is now being solved on a finer grid */ 889 PetscCall(PCMGSetLevels(pc,levels,NULL)); 890 n = levels; 891 PetscCall(PCSetFromOptions(pc)); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 892 mglevels = mg->levels; 893 } 894 } 895 PetscCall(KSPGetPC(mglevels[0]->smoothd,&cpc)); 896 897 /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 898 /* so use those from global PC */ 899 /* Is this what we always want? What if user wants to keep old one? */ 900 PetscCall(KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset)); 901 if (opsset) { 902 Mat mmat; 903 PetscCall(KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat)); 904 if (mmat == pc->pmat) opsset = PETSC_FALSE; 905 } 906 907 /* Create CR solvers */ 908 PetscCall(PCMGGetAdaptCR(pc, &doCR)); 909 if (doCR) { 910 const char *prefix; 911 912 PetscCall(PCGetOptionsPrefix(pc, &prefix)); 913 for (i = 1; i < n; ++i) { 914 PC ipc, cr; 915 char crprefix[128]; 916 917 PetscCall(KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr)); 918 PetscCall(KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE)); 919 PetscCall(PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i)); 920 PetscCall(KSPSetOptionsPrefix(mglevels[i]->cr, prefix)); 921 PetscCall(PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level)); 922 PetscCall(KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV)); 923 PetscCall(KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL)); 924 PetscCall(KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED)); 925 PetscCall(KSPGetPC(mglevels[i]->cr, &ipc)); 926 927 PetscCall(PCSetType(ipc, PCCOMPOSITE)); 928 PetscCall(PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE)); 929 PetscCall(PCCompositeAddPCType(ipc, PCSOR)); 930 PetscCall(CreateCR_Private(pc, i, &cr)); 931 PetscCall(PCCompositeAddPC(ipc, cr)); 932 PetscCall(PCDestroy(&cr)); 933 934 PetscCall(KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd)); 935 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE)); 936 PetscCall(PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i)); 937 PetscCall(KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix)); 938 PetscCall(PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr)); 939 } 940 } 941 942 if (!opsset) { 943 PetscCall(PCGetUseAmat(pc,&use_amat)); 944 if (use_amat) { 945 PetscCall(PetscInfo(pc,"Using outer operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n")); 946 PetscCall(KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat)); 947 } else { 948 PetscCall(PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n")); 949 PetscCall(KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat)); 950 } 951 } 952 953 for (i=n-1; i>0; i--) { 954 if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 955 missinginterpolate = PETSC_TRUE; 956 break; 957 } 958 } 959 960 PetscCall(KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB)); 961 if (dA == dB) dAeqdB = PETSC_TRUE; 962 if (mg->galerkin == PC_MG_GALERKIN_NONE || ((mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_MAT) && !dAeqdB)) { 963 needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 964 } 965 966 if (pc->dm && !pc->setupcalled) { 967 /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 968 PetscCall(KSPSetDM(mglevels[n-1]->smoothd,pc->dm)); 969 PetscCall(KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE)); 970 if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) { 971 PetscCall(KSPSetDM(mglevels[n-1]->smoothu,pc->dm)); 972 PetscCall(KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE)); 973 } 974 if (mglevels[n-1]->cr) { 975 PetscCall(KSPSetDM(mglevels[n-1]->cr,pc->dm)); 976 PetscCall(KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE)); 977 } 978 } 979 980 /* 981 Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 982 Skipping for externally managed hierarchy (such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 983 */ 984 if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 985 /* construct the interpolation from the DMs */ 986 Mat p; 987 Vec rscale; 988 PetscCall(PetscMalloc1(n,&dms)); 989 dms[n-1] = pc->dm; 990 /* Separately create them so we do not get DMKSP interference between levels */ 991 for (i=n-2; i>-1; i--) PetscCall(DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i])); 992 for (i=n-2; i>-1; i--) { 993 DMKSP kdm; 994 PetscBool dmhasrestrict, dmhasinject; 995 PetscCall(KSPSetDM(mglevels[i]->smoothd,dms[i])); 996 if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE)); 997 if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 998 PetscCall(KSPSetDM(mglevels[i]->smoothu,dms[i])); 999 if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE)); 1000 } 1001 if (mglevels[i]->cr) { 1002 PetscCall(KSPSetDM(mglevels[i]->cr,dms[i])); 1003 if (!needRestricts) PetscCall(KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE)); 1004 } 1005 PetscCall(DMGetDMKSPWrite(dms[i],&kdm)); 1006 /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 1007 * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 1008 kdm->ops->computerhs = NULL; 1009 kdm->rhsctx = NULL; 1010 if (!mglevels[i+1]->interpolate) { 1011 PetscCall(DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale)); 1012 PetscCall(PCMGSetInterpolation(pc,i+1,p)); 1013 if (rscale) PetscCall(PCMGSetRScale(pc,i+1,rscale)); 1014 PetscCall(VecDestroy(&rscale)); 1015 PetscCall(MatDestroy(&p)); 1016 } 1017 PetscCall(DMHasCreateRestriction(dms[i],&dmhasrestrict)); 1018 if (dmhasrestrict && !mglevels[i+1]->restrct) { 1019 PetscCall(DMCreateRestriction(dms[i],dms[i+1],&p)); 1020 PetscCall(PCMGSetRestriction(pc,i+1,p)); 1021 PetscCall(MatDestroy(&p)); 1022 } 1023 PetscCall(DMHasCreateInjection(dms[i],&dmhasinject)); 1024 if (dmhasinject && !mglevels[i+1]->inject) { 1025 PetscCall(DMCreateInjection(dms[i],dms[i+1],&p)); 1026 PetscCall(PCMGSetInjection(pc,i+1,p)); 1027 PetscCall(MatDestroy(&p)); 1028 } 1029 } 1030 1031 for (i=n-2; i>-1; i--) PetscCall(DMDestroy(&dms[i])); 1032 PetscCall(PetscFree(dms)); 1033 } 1034 1035 if (mg->galerkin < PC_MG_GALERKIN_NONE) { 1036 Mat A,B; 1037 PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 1038 MatReuse reuse = MAT_INITIAL_MATRIX; 1039 1040 if (mg->galerkin == PC_MG_GALERKIN_PMAT || mg->galerkin == PC_MG_GALERKIN_BOTH) doB = PETSC_TRUE; 1041 if (mg->galerkin == PC_MG_GALERKIN_MAT || (mg->galerkin == PC_MG_GALERKIN_BOTH && dA != dB)) doA = PETSC_TRUE; 1042 if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 1043 for (i=n-2; i>-1; i--) { 1044 if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) { /* see if we can compute a coarse space */ 1045 PetscCall(PCMGComputeCoarseSpace_Internal(pc, i+1, mg->coarseSpaceType, mg->Nc, NULL, &mglevels[i+1]->coarseSpace)); 1046 PetscCall(PCMGSetInterpolation(pc,i+1,mglevels[i+1]->coarseSpace)); 1047 } 1048 PetscCheck(mglevels[i+1]->restrct || mglevels[i+1]->interpolate,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0"); 1049 if (!mglevels[i+1]->interpolate) { 1050 PetscCall(PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct)); 1051 } 1052 if (!mglevels[i+1]->restrct) { 1053 PetscCall(PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate)); 1054 } 1055 if (reuse == MAT_REUSE_MATRIX) { 1056 PetscCall(KSPGetOperators(mglevels[i]->smoothd,&A,&B)); 1057 } 1058 if (doA) { 1059 PetscCall(MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A)); 1060 } 1061 if (doB) { 1062 PetscCall(MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B)); 1063 } 1064 /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 1065 if (!doA && dAeqdB) { 1066 if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)B)); 1067 A = B; 1068 } else if (!doA && reuse == MAT_INITIAL_MATRIX) { 1069 PetscCall(KSPGetOperators(mglevels[i]->smoothd,&A,NULL)); 1070 PetscCall(PetscObjectReference((PetscObject)A)); 1071 } 1072 if (!doB && dAeqdB) { 1073 if (reuse == MAT_INITIAL_MATRIX) PetscCall(PetscObjectReference((PetscObject)A)); 1074 B = A; 1075 } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 1076 PetscCall(KSPGetOperators(mglevels[i]->smoothd,NULL,&B)); 1077 PetscCall(PetscObjectReference((PetscObject)B)); 1078 } 1079 if (reuse == MAT_INITIAL_MATRIX) { 1080 PetscCall(KSPSetOperators(mglevels[i]->smoothd,A,B)); 1081 PetscCall(PetscObjectDereference((PetscObject)A)); 1082 PetscCall(PetscObjectDereference((PetscObject)B)); 1083 } 1084 dA = A; 1085 dB = B; 1086 } 1087 } 1088 1089 /* Adapt interpolation matrices */ 1090 if (adaptInterpolation) { 1091 for (i = 0; i < n; ++i) { 1092 if (!mglevels[i]->coarseSpace) { 1093 PetscCall(PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace)); 1094 } 1095 if (i) PetscCall(PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace)); 1096 } 1097 for (i = n-2; i > -1; --i) { 1098 PetscCall(PCMGRecomputeLevelOperators_Internal(pc, i)); 1099 } 1100 } 1101 1102 if (needRestricts && pc->dm) { 1103 for (i=n-2; i>=0; i--) { 1104 DM dmfine,dmcoarse; 1105 Mat Restrict,Inject; 1106 Vec rscale; 1107 PetscCall(KSPGetDM(mglevels[i+1]->smoothd,&dmfine)); 1108 PetscCall(KSPGetDM(mglevels[i]->smoothd,&dmcoarse)); 1109 PetscCall(PCMGGetRestriction(pc,i+1,&Restrict)); 1110 PetscCall(PCMGGetRScale(pc,i+1,&rscale)); 1111 PetscCall(PCMGGetInjection(pc,i+1,&Inject)); 1112 PetscCall(DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse)); 1113 } 1114 } 1115 1116 if (!pc->setupcalled) { 1117 for (i=0; i<n; i++) { 1118 PetscCall(KSPSetFromOptions(mglevels[i]->smoothd)); 1119 } 1120 for (i=1; i<n; i++) { 1121 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 1122 PetscCall(KSPSetFromOptions(mglevels[i]->smoothu)); 1123 } 1124 if (mglevels[i]->cr) { 1125 PetscCall(KSPSetFromOptions(mglevels[i]->cr)); 1126 } 1127 } 1128 /* insure that if either interpolation or restriction is set the other other one is set */ 1129 for (i=1; i<n; i++) { 1130 PetscCall(PCMGGetInterpolation(pc,i,NULL)); 1131 PetscCall(PCMGGetRestriction(pc,i,NULL)); 1132 } 1133 for (i=0; i<n-1; i++) { 1134 if (!mglevels[i]->b) { 1135 Vec *vec; 1136 PetscCall(KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL)); 1137 PetscCall(PCMGSetRhs(pc,i,*vec)); 1138 PetscCall(VecDestroy(vec)); 1139 PetscCall(PetscFree(vec)); 1140 } 1141 if (!mglevels[i]->r && i) { 1142 PetscCall(VecDuplicate(mglevels[i]->b,&tvec)); 1143 PetscCall(PCMGSetR(pc,i,tvec)); 1144 PetscCall(VecDestroy(&tvec)); 1145 } 1146 if (!mglevels[i]->x) { 1147 PetscCall(VecDuplicate(mglevels[i]->b,&tvec)); 1148 PetscCall(PCMGSetX(pc,i,tvec)); 1149 PetscCall(VecDestroy(&tvec)); 1150 } 1151 if (doCR) { 1152 PetscCall(VecDuplicate(mglevels[i]->b,&mglevels[i]->crx)); 1153 PetscCall(VecDuplicate(mglevels[i]->b,&mglevels[i]->crb)); 1154 PetscCall(VecZeroEntries(mglevels[i]->crb)); 1155 } 1156 } 1157 if (n != 1 && !mglevels[n-1]->r) { 1158 /* PCMGSetR() on the finest level if user did not supply it */ 1159 Vec *vec; 1160 PetscCall(KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL)); 1161 PetscCall(PCMGSetR(pc,n-1,*vec)); 1162 PetscCall(VecDestroy(vec)); 1163 PetscCall(PetscFree(vec)); 1164 } 1165 if (doCR) { 1166 PetscCall(VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx)); 1167 PetscCall(VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb)); 1168 PetscCall(VecZeroEntries(mglevels[n-1]->crb)); 1169 } 1170 } 1171 1172 if (pc->dm) { 1173 /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 1174 for (i=0; i<n-1; i++) { 1175 if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 1176 } 1177 } 1178 // We got here (PCSetUp_MG) because the matrix has changed, which means the smoother needs to be set up again (e.g., 1179 // new diagonal for Jacobi). Setting it here allows it to be logged under PCSetUp rather than deep inside a PCApply. 1180 if (mglevels[n-1]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[n-1]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 1181 1182 for (i=1; i<n; i++) { 1183 if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1) { 1184 /* if doing only down then initial guess is zero */ 1185 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE)); 1186 } 1187 if (mglevels[i]->cr) PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE)); 1188 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0)); 1189 PetscCall(KSPSetUp(mglevels[i]->smoothd)); 1190 if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1191 pc->failedreason = PC_SUBPC_ERROR; 1192 } 1193 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0)); 1194 if (!mglevels[i]->residual) { 1195 Mat mat; 1196 PetscCall(KSPGetOperators(mglevels[i]->smoothd,&mat,NULL)); 1197 PetscCall(PCMGSetResidual(pc,i,PCMGResidualDefault,mat)); 1198 } 1199 if (!mglevels[i]->residualtranspose) { 1200 Mat mat; 1201 PetscCall(KSPGetOperators(mglevels[i]->smoothd,&mat,NULL)); 1202 PetscCall(PCMGSetResidualTranspose(pc,i,PCMGResidualTransposeDefault,mat)); 1203 } 1204 } 1205 for (i=1; i<n; i++) { 1206 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 1207 Mat downmat,downpmat; 1208 1209 /* check if operators have been set for up, if not use down operators to set them */ 1210 PetscCall(KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL)); 1211 if (!opsset) { 1212 PetscCall(KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat)); 1213 PetscCall(KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat)); 1214 } 1215 1216 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE)); 1217 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0)); 1218 PetscCall(KSPSetUp(mglevels[i]->smoothu)); 1219 if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) { 1220 pc->failedreason = PC_SUBPC_ERROR; 1221 } 1222 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0)); 1223 } 1224 if (mglevels[i]->cr) { 1225 Mat downmat,downpmat; 1226 1227 /* check if operators have been set for up, if not use down operators to set them */ 1228 PetscCall(KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL)); 1229 if (!opsset) { 1230 PetscCall(KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat)); 1231 PetscCall(KSPSetOperators(mglevels[i]->cr,downmat,downpmat)); 1232 } 1233 1234 PetscCall(KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE)); 1235 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0)); 1236 PetscCall(KSPSetUp(mglevels[i]->cr)); 1237 if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) { 1238 pc->failedreason = PC_SUBPC_ERROR; 1239 } 1240 if (mglevels[i]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0)); 1241 } 1242 } 1243 1244 if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0)); 1245 PetscCall(KSPSetUp(mglevels[0]->smoothd)); 1246 if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1247 pc->failedreason = PC_SUBPC_ERROR; 1248 } 1249 if (mglevels[0]->eventsmoothsetup) PetscCall(PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0)); 1250 1251 /* 1252 Dump the interpolation/restriction matrices plus the 1253 Jacobian/stiffness on each level. This allows MATLAB users to 1254 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 1255 1256 Only support one or the other at the same time. 1257 */ 1258 #if defined(PETSC_USE_SOCKET_VIEWER) 1259 PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL)); 1260 if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 1261 dump = PETSC_FALSE; 1262 #endif 1263 PetscCall(PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL)); 1264 if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 1265 1266 if (viewer) { 1267 for (i=1; i<n; i++) { 1268 PetscCall(MatView(mglevels[i]->restrct,viewer)); 1269 } 1270 for (i=0; i<n; i++) { 1271 PetscCall(KSPGetPC(mglevels[i]->smoothd,&pc)); 1272 PetscCall(MatView(pc->mat,viewer)); 1273 } 1274 } 1275 PetscFunctionReturn(0); 1276 } 1277 1278 /* -------------------------------------------------------------------------------------*/ 1279 1280 PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 1281 { 1282 PC_MG *mg = (PC_MG *) pc->data; 1283 1284 PetscFunctionBegin; 1285 *levels = mg->nlevels; 1286 PetscFunctionReturn(0); 1287 } 1288 1289 /*@ 1290 PCMGGetLevels - Gets the number of levels to use with MG. 1291 1292 Not Collective 1293 1294 Input Parameter: 1295 . pc - the preconditioner context 1296 1297 Output parameter: 1298 . levels - the number of levels 1299 1300 Level: advanced 1301 1302 .seealso: `PCMGSetLevels()` 1303 @*/ 1304 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 1305 { 1306 PetscFunctionBegin; 1307 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1308 PetscValidIntPointer(levels,2); 1309 *levels = 0; 1310 PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels)); 1311 PetscFunctionReturn(0); 1312 } 1313 1314 /*@ 1315 PCMGGetGridComplexity - compute operator and grid complexity of MG hierarchy 1316 1317 Input Parameter: 1318 . pc - the preconditioner context 1319 1320 Output Parameters: 1321 + gc - grid complexity = sum_i(n_i) / n_0 1322 - oc - operator complexity = sum_i(nnz_i) / nnz_0 1323 1324 Level: advanced 1325 1326 .seealso: `PCMGGetLevels()` 1327 @*/ 1328 PetscErrorCode PCMGGetGridComplexity(PC pc, PetscReal *gc, PetscReal *oc) 1329 { 1330 PC_MG *mg = (PC_MG*)pc->data; 1331 PC_MG_Levels **mglevels = mg->levels; 1332 PetscInt lev,N; 1333 PetscLogDouble nnz0 = 0, sgc = 0, soc = 0, n0 = 0; 1334 MatInfo info; 1335 1336 PetscFunctionBegin; 1337 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1338 if (gc) PetscValidRealPointer(gc,2); 1339 if (oc) PetscValidRealPointer(oc,3); 1340 if (!pc->setupcalled) { 1341 if (gc) *gc = 0; 1342 if (oc) *oc = 0; 1343 PetscFunctionReturn(0); 1344 } 1345 PetscCheck(mg->nlevels > 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MG has no levels"); 1346 for (lev=0; lev<mg->nlevels; lev++) { 1347 Mat dB; 1348 PetscCall(KSPGetOperators(mglevels[lev]->smoothd,NULL,&dB)); 1349 PetscCall(MatGetInfo(dB,MAT_GLOBAL_SUM,&info)); /* global reduction */ 1350 PetscCall(MatGetSize(dB,&N,NULL)); 1351 sgc += N; 1352 soc += info.nz_used; 1353 if (lev==mg->nlevels-1) {nnz0 = info.nz_used; n0 = N;} 1354 } 1355 if (n0 > 0 && gc) *gc = (PetscReal)(sgc/n0); 1356 else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number for grid points on finest level is not available"); 1357 if (nnz0 > 0 && oc) *oc = (PetscReal)(soc/nnz0); 1358 PetscFunctionReturn(0); 1359 } 1360 1361 /*@ 1362 PCMGSetType - Determines the form of multigrid to use: 1363 multiplicative, additive, full, or the Kaskade algorithm. 1364 1365 Logically Collective on PC 1366 1367 Input Parameters: 1368 + pc - the preconditioner context 1369 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 1370 PC_MG_FULL, PC_MG_KASKADE 1371 1372 Options Database Key: 1373 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 1374 additive, full, kaskade 1375 1376 Level: advanced 1377 1378 .seealso: `PCMGSetLevels()` 1379 @*/ 1380 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 1381 { 1382 PC_MG *mg = (PC_MG*)pc->data; 1383 1384 PetscFunctionBegin; 1385 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1386 PetscValidLogicalCollectiveEnum(pc,form,2); 1387 mg->am = form; 1388 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 1389 else pc->ops->applyrichardson = NULL; 1390 PetscFunctionReturn(0); 1391 } 1392 1393 /*@ 1394 PCMGGetType - Determines the form of multigrid to use: 1395 multiplicative, additive, full, or the Kaskade algorithm. 1396 1397 Logically Collective on PC 1398 1399 Input Parameter: 1400 . pc - the preconditioner context 1401 1402 Output Parameter: 1403 . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 1404 1405 Level: advanced 1406 1407 .seealso: `PCMGSetLevels()` 1408 @*/ 1409 PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 1410 { 1411 PC_MG *mg = (PC_MG*)pc->data; 1412 1413 PetscFunctionBegin; 1414 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1415 *type = mg->am; 1416 PetscFunctionReturn(0); 1417 } 1418 1419 /*@ 1420 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 1421 complicated cycling. 1422 1423 Logically Collective on PC 1424 1425 Input Parameters: 1426 + pc - the multigrid context 1427 - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 1428 1429 Options Database Key: 1430 . -pc_mg_cycle_type <v,w> - provide the cycle desired 1431 1432 Level: advanced 1433 1434 .seealso: `PCMGSetCycleTypeOnLevel()` 1435 @*/ 1436 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 1437 { 1438 PC_MG *mg = (PC_MG*)pc->data; 1439 PC_MG_Levels **mglevels = mg->levels; 1440 PetscInt i,levels; 1441 1442 PetscFunctionBegin; 1443 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1444 PetscValidLogicalCollectiveEnum(pc,n,2); 1445 PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1446 levels = mglevels[0]->levels; 1447 for (i=0; i<levels; i++) mglevels[i]->cycles = n; 1448 PetscFunctionReturn(0); 1449 } 1450 1451 /*@ 1452 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 1453 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 1454 1455 Logically Collective on PC 1456 1457 Input Parameters: 1458 + pc - the multigrid context 1459 - n - number of cycles (default is 1) 1460 1461 Options Database Key: 1462 . -pc_mg_multiplicative_cycles n - set the number of cycles 1463 1464 Level: advanced 1465 1466 Notes: 1467 This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 1468 1469 .seealso: `PCMGSetCycleTypeOnLevel()`, `PCMGSetCycleType()` 1470 @*/ 1471 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 1472 { 1473 PC_MG *mg = (PC_MG*)pc->data; 1474 1475 PetscFunctionBegin; 1476 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1477 PetscValidLogicalCollectiveInt(pc,n,2); 1478 mg->cyclesperpcapply = n; 1479 PetscFunctionReturn(0); 1480 } 1481 1482 PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1483 { 1484 PC_MG *mg = (PC_MG*)pc->data; 1485 1486 PetscFunctionBegin; 1487 mg->galerkin = use; 1488 PetscFunctionReturn(0); 1489 } 1490 1491 /*@ 1492 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 1493 finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1494 1495 Logically Collective on PC 1496 1497 Input Parameters: 1498 + pc - the multigrid context 1499 - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1500 1501 Options Database Key: 1502 . -pc_mg_galerkin <both,pmat,mat,none> - set the matrices to form via the Galerkin process 1503 1504 Level: intermediate 1505 1506 Notes: 1507 Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 1508 use the PCMG construction of the coarser grids. 1509 1510 .seealso: `PCMGGetGalerkin()`, `PCMGGalerkinType` 1511 1512 @*/ 1513 PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1514 { 1515 PetscFunctionBegin; 1516 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1517 PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use)); 1518 PetscFunctionReturn(0); 1519 } 1520 1521 /*@ 1522 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 1523 A_i-1 = r_i * A_i * p_i 1524 1525 Not Collective 1526 1527 Input Parameter: 1528 . pc - the multigrid context 1529 1530 Output Parameter: 1531 . galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL 1532 1533 Level: intermediate 1534 1535 .seealso: `PCMGSetGalerkin()`, `PCMGGalerkinType` 1536 1537 @*/ 1538 PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 1539 { 1540 PC_MG *mg = (PC_MG*)pc->data; 1541 1542 PetscFunctionBegin; 1543 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1544 *galerkin = mg->galerkin; 1545 PetscFunctionReturn(0); 1546 } 1547 1548 PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt) 1549 { 1550 PC_MG *mg = (PC_MG *) pc->data; 1551 1552 PetscFunctionBegin; 1553 mg->adaptInterpolation = adapt; 1554 PetscFunctionReturn(0); 1555 } 1556 1557 PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt) 1558 { 1559 PC_MG *mg = (PC_MG *) pc->data; 1560 1561 PetscFunctionBegin; 1562 *adapt = mg->adaptInterpolation; 1563 PetscFunctionReturn(0); 1564 } 1565 1566 PetscErrorCode PCMGSetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType ctype) 1567 { 1568 PC_MG *mg = (PC_MG *) pc->data; 1569 1570 PetscFunctionBegin; 1571 mg->adaptInterpolation = ctype != PCMG_ADAPT_NONE ? PETSC_TRUE : PETSC_FALSE; 1572 mg->coarseSpaceType = ctype; 1573 PetscFunctionReturn(0); 1574 } 1575 1576 PetscErrorCode PCMGGetAdaptCoarseSpaceType_MG(PC pc, PCMGCoarseSpaceType *ctype) 1577 { 1578 PC_MG *mg = (PC_MG *) pc->data; 1579 1580 PetscFunctionBegin; 1581 *ctype = mg->coarseSpaceType; 1582 PetscFunctionReturn(0); 1583 } 1584 1585 PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr) 1586 { 1587 PC_MG *mg = (PC_MG *) pc->data; 1588 1589 PetscFunctionBegin; 1590 mg->compatibleRelaxation = cr; 1591 PetscFunctionReturn(0); 1592 } 1593 1594 PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr) 1595 { 1596 PC_MG *mg = (PC_MG *) pc->data; 1597 1598 PetscFunctionBegin; 1599 *cr = mg->compatibleRelaxation; 1600 PetscFunctionReturn(0); 1601 } 1602 1603 /*@C 1604 PCMGSetAdaptCoarseSpaceType - Set the type of adaptive coarse space. 1605 1606 Adapts or creates the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1607 1608 Logically Collective on PC 1609 1610 Input Parameters: 1611 + pc - the multigrid context 1612 - ctype - the type of coarse space 1613 1614 Options Database Keys: 1615 + -pc_mg_adapt_interp_n <int> - The number of modes to use 1616 - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: none, polynomial, harmonic, eigenvector, generalized_eigenvector, gdsw 1617 1618 Level: intermediate 1619 1620 .keywords: MG, set, Galerkin 1621 .seealso: `PCMGCoarseSpaceType`, `PCMGGetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()` 1622 @*/ 1623 PetscErrorCode PCMGSetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType ctype) 1624 { 1625 PetscFunctionBegin; 1626 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1627 PetscValidLogicalCollectiveEnum(pc,ctype,2); 1628 PetscTryMethod(pc,"PCMGSetAdaptCoarseSpaceType_C",(PC,PCMGCoarseSpaceType),(pc,ctype)); 1629 PetscFunctionReturn(0); 1630 } 1631 1632 /*@C 1633 PCMGGetAdaptCoarseSpaceType - Get the type of adaptive coarse space. 1634 1635 Not Collective 1636 1637 Input Parameter: 1638 . pc - the multigrid context 1639 1640 Output Parameter: 1641 . ctype - the type of coarse space 1642 1643 Level: intermediate 1644 1645 .keywords: MG, Get, Galerkin 1646 .seealso: `PCMGCoarseSpaceType`, `PCMGSetAdaptCoarseSpaceType()`, `PCMGSetGalerkin()` 1647 @*/ 1648 PetscErrorCode PCMGGetAdaptCoarseSpaceType(PC pc, PCMGCoarseSpaceType *ctype) 1649 { 1650 PetscFunctionBegin; 1651 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1652 PetscValidPointer(ctype,2); 1653 PetscUseMethod(pc,"PCMGGetAdaptCoarseSpaceType_C",(PC,PCMGCoarseSpaceType*),(pc,ctype)); 1654 PetscFunctionReturn(0); 1655 } 1656 1657 /* MATT: REMOVE? */ 1658 /*@ 1659 PCMGSetAdaptInterpolation - Adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1660 1661 Logically Collective on PC 1662 1663 Input Parameters: 1664 + pc - the multigrid context 1665 - adapt - flag for adaptation of the interpolator 1666 1667 Options Database Keys: 1668 + -pc_mg_adapt_interp - Turn on adaptation 1669 . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension 1670 - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector 1671 1672 Level: intermediate 1673 1674 .keywords: MG, set, Galerkin 1675 .seealso: `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()` 1676 @*/ 1677 PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt) 1678 { 1679 PetscFunctionBegin; 1680 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1681 PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt)); 1682 PetscFunctionReturn(0); 1683 } 1684 1685 /*@ 1686 PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1687 1688 Not collective 1689 1690 Input Parameter: 1691 . pc - the multigrid context 1692 1693 Output Parameter: 1694 . adapt - flag for adaptation of the interpolator 1695 1696 Level: intermediate 1697 1698 .keywords: MG, set, Galerkin 1699 .seealso: `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()` 1700 @*/ 1701 PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt) 1702 { 1703 PetscFunctionBegin; 1704 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1705 PetscValidBoolPointer(adapt, 2); 1706 PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt)); 1707 PetscFunctionReturn(0); 1708 } 1709 1710 /*@ 1711 PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation. 1712 1713 Logically Collective on PC 1714 1715 Input Parameters: 1716 + pc - the multigrid context 1717 - cr - flag for compatible relaxation 1718 1719 Options Database Keys: 1720 . -pc_mg_adapt_cr - Turn on compatible relaxation 1721 1722 Level: intermediate 1723 1724 .keywords: MG, set, Galerkin 1725 .seealso: `PCMGGetAdaptCR()`, `PCMGSetAdaptInterpolation()`, `PCMGSetGalerkin()` 1726 @*/ 1727 PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr) 1728 { 1729 PetscFunctionBegin; 1730 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1731 PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr)); 1732 PetscFunctionReturn(0); 1733 } 1734 1735 /*@ 1736 PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation. 1737 1738 Not collective 1739 1740 Input Parameter: 1741 . pc - the multigrid context 1742 1743 Output Parameter: 1744 . cr - flag for compatible relaxaion 1745 1746 Level: intermediate 1747 1748 .keywords: MG, set, Galerkin 1749 .seealso: `PCMGSetAdaptCR()`, `PCMGGetAdaptInterpolation()`, `PCMGSetGalerkin()` 1750 @*/ 1751 PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr) 1752 { 1753 PetscFunctionBegin; 1754 PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1755 PetscValidBoolPointer(cr, 2); 1756 PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr)); 1757 PetscFunctionReturn(0); 1758 } 1759 1760 /*@ 1761 PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1762 on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1763 pre- and post-smoothing steps. 1764 1765 Logically Collective on PC 1766 1767 Input Parameters: 1768 + mg - the multigrid context 1769 - n - the number of smoothing steps 1770 1771 Options Database Key: 1772 . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 1773 1774 Level: advanced 1775 1776 Notes: 1777 this does not set a value on the coarsest grid, since we assume that 1778 there is no separate smooth up on the coarsest grid. 1779 1780 .seealso: `PCMGSetDistinctSmoothUp()` 1781 @*/ 1782 PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 1783 { 1784 PC_MG *mg = (PC_MG*)pc->data; 1785 PC_MG_Levels **mglevels = mg->levels; 1786 PetscInt i,levels; 1787 1788 PetscFunctionBegin; 1789 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1790 PetscValidLogicalCollectiveInt(pc,n,2); 1791 PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1792 levels = mglevels[0]->levels; 1793 1794 for (i=1; i<levels; i++) { 1795 PetscCall(KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n)); 1796 PetscCall(KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n)); 1797 mg->default_smoothu = n; 1798 mg->default_smoothd = n; 1799 } 1800 PetscFunctionReturn(0); 1801 } 1802 1803 /*@ 1804 PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels 1805 and adds the suffix _up to the options name 1806 1807 Logically Collective on PC 1808 1809 Input Parameters: 1810 . pc - the preconditioner context 1811 1812 Options Database Key: 1813 . -pc_mg_distinct_smoothup <bool> - use distinct smoothing objects 1814 1815 Level: advanced 1816 1817 Notes: 1818 this does not set a value on the coarsest grid, since we assume that 1819 there is no separate smooth up on the coarsest grid. 1820 1821 .seealso: `PCMGSetNumberSmooth()` 1822 @*/ 1823 PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1824 { 1825 PC_MG *mg = (PC_MG*)pc->data; 1826 PC_MG_Levels **mglevels = mg->levels; 1827 PetscInt i,levels; 1828 KSP subksp; 1829 1830 PetscFunctionBegin; 1831 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1832 PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1833 levels = mglevels[0]->levels; 1834 1835 for (i=1; i<levels; i++) { 1836 const char *prefix = NULL; 1837 /* make sure smoother up and down are different */ 1838 PetscCall(PCMGGetSmootherUp(pc,i,&subksp)); 1839 PetscCall(KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix)); 1840 PetscCall(KSPSetOptionsPrefix(subksp,prefix)); 1841 PetscCall(KSPAppendOptionsPrefix(subksp,"up_")); 1842 } 1843 PetscFunctionReturn(0); 1844 } 1845 1846 /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 1847 PetscErrorCode PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[]) 1848 { 1849 PC_MG *mg = (PC_MG*)pc->data; 1850 PC_MG_Levels **mglevels = mg->levels; 1851 Mat *mat; 1852 PetscInt l; 1853 1854 PetscFunctionBegin; 1855 PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1856 PetscCall(PetscMalloc1(mg->nlevels,&mat)); 1857 for (l=1; l< mg->nlevels; l++) { 1858 mat[l-1] = mglevels[l]->interpolate; 1859 PetscCall(PetscObjectReference((PetscObject)mat[l-1])); 1860 } 1861 *num_levels = mg->nlevels; 1862 *interpolations = mat; 1863 PetscFunctionReturn(0); 1864 } 1865 1866 /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 1867 PetscErrorCode PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[]) 1868 { 1869 PC_MG *mg = (PC_MG*)pc->data; 1870 PC_MG_Levels **mglevels = mg->levels; 1871 PetscInt l; 1872 Mat *mat; 1873 1874 PetscFunctionBegin; 1875 PetscCheck(mglevels,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 1876 PetscCall(PetscMalloc1(mg->nlevels,&mat)); 1877 for (l=0; l<mg->nlevels-1; l++) { 1878 PetscCall(KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]))); 1879 PetscCall(PetscObjectReference((PetscObject)mat[l])); 1880 } 1881 *num_levels = mg->nlevels; 1882 *coarseOperators = mat; 1883 PetscFunctionReturn(0); 1884 } 1885 1886 /*@C 1887 PCMGRegisterCoarseSpaceConstructor - Adds a method to the PCMG package for coarse space construction. 1888 1889 Not collective 1890 1891 Input Parameters: 1892 + name - name of the constructor 1893 - function - constructor routine 1894 1895 Notes: 1896 Calling sequence for the routine: 1897 $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp) 1898 $ pc - The PC object 1899 $ l - The multigrid level, 0 is the coarse level 1900 $ dm - The DM for this level 1901 $ smooth - The level smoother 1902 $ Nc - The size of the coarse space 1903 $ initGuess - Basis for an initial guess for the space 1904 $ coarseSp - A basis for the computed coarse space 1905 1906 Level: advanced 1907 1908 .seealso: `PCMGGetCoarseSpaceConstructor()`, `PCRegister()` 1909 @*/ 1910 PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*)) 1911 { 1912 PetscFunctionBegin; 1913 PetscCall(PCInitializePackage()); 1914 PetscCall(PetscFunctionListAdd(&PCMGCoarseList,name,function)); 1915 PetscFunctionReturn(0); 1916 } 1917 1918 /*@C 1919 PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method. 1920 1921 Not collective 1922 1923 Input Parameter: 1924 . name - name of the constructor 1925 1926 Output Parameter: 1927 . function - constructor routine 1928 1929 Notes: 1930 Calling sequence for the routine: 1931 $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp) 1932 $ pc - The PC object 1933 $ l - The multigrid level, 0 is the coarse level 1934 $ dm - The DM for this level 1935 $ smooth - The level smoother 1936 $ Nc - The size of the coarse space 1937 $ initGuess - Basis for an initial guess for the space 1938 $ coarseSp - A basis for the computed coarse space 1939 1940 Level: advanced 1941 1942 .seealso: `PCMGRegisterCoarseSpaceConstructor()`, `PCRegister()` 1943 @*/ 1944 PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*)) 1945 { 1946 PetscFunctionBegin; 1947 PetscCall(PetscFunctionListFind(PCMGCoarseList,name,function)); 1948 PetscFunctionReturn(0); 1949 } 1950 1951 /* ----------------------------------------------------------------------------------------*/ 1952 1953 /*MC 1954 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1955 information about the coarser grid matrices and restriction/interpolation operators. 1956 1957 Options Database Keys: 1958 + -pc_mg_levels <nlevels> - number of levels including finest 1959 . -pc_mg_cycle_type <v,w> - provide the cycle desired 1960 . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 1961 . -pc_mg_log - log information about time spent on each level of the solver 1962 . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 1963 . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 1964 . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 1965 . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1966 to the Socket viewer for reading from MATLAB. 1967 - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 1968 to the binary output file called binaryoutput 1969 1970 Notes: 1971 If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method 1972 1973 When run with a single level the smoother options are used on that level NOT the coarse grid solver options 1974 1975 When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 1976 is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 1977 (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 1978 residual is computed at the end of each cycle. 1979 1980 Level: intermediate 1981 1982 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCMGType`, `PCEXOTIC`, `PCGAMG`, `PCML`, `PCHYPRE` 1983 `PCMGSetLevels()`, `PCMGGetLevels()`, `PCMGSetType()`, `PCMGSetCycleType()`, 1984 `PCMGSetDistinctSmoothUp()`, `PCMGGetCoarseSolve()`, `PCMGSetResidual()`, `PCMGSetInterpolation()`, 1985 `PCMGSetRestriction()`, `PCMGGetSmoother()`, `PCMGGetSmootherUp()`, `PCMGGetSmootherDown()`, 1986 `PCMGSetCycleTypeOnLevel()`, `PCMGSetRhs()`, `PCMGSetX()`, `PCMGSetR()` 1987 M*/ 1988 1989 PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 1990 { 1991 PC_MG *mg; 1992 1993 PetscFunctionBegin; 1994 PetscCall(PetscNewLog(pc,&mg)); 1995 pc->data = mg; 1996 mg->nlevels = -1; 1997 mg->am = PC_MG_MULTIPLICATIVE; 1998 mg->galerkin = PC_MG_GALERKIN_NONE; 1999 mg->adaptInterpolation = PETSC_FALSE; 2000 mg->Nc = -1; 2001 mg->eigenvalue = -1; 2002 2003 pc->useAmat = PETSC_TRUE; 2004 2005 pc->ops->apply = PCApply_MG; 2006 pc->ops->applytranspose = PCApplyTranspose_MG; 2007 pc->ops->matapply = PCMatApply_MG; 2008 pc->ops->setup = PCSetUp_MG; 2009 pc->ops->reset = PCReset_MG; 2010 pc->ops->destroy = PCDestroy_MG; 2011 pc->ops->setfromoptions = PCSetFromOptions_MG; 2012 pc->ops->view = PCView_MG; 2013 2014 PetscCall(PetscObjectComposedDataRegister(&mg->eigenvalue)); 2015 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG)); 2016 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG)); 2017 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG)); 2018 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG)); 2019 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG)); 2020 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG)); 2021 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG)); 2022 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG)); 2023 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG)); 2024 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCoarseSpaceType_C",PCMGSetAdaptCoarseSpaceType_MG)); 2025 PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCoarseSpaceType_C",PCMGGetAdaptCoarseSpaceType_MG)); 2026 PetscFunctionReturn(0); 2027 } 2028