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