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