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