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