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 } 634 else { 635 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); 636 ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 637 } 638 } 639 640 for (i=n-1; i>0; i--) { 641 if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 642 missinginterpolate = PETSC_TRUE; 643 continue; 644 } 645 } 646 647 ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 648 if (dA == dB) dAeqdB = PETSC_TRUE; 649 if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) { 650 needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 651 } 652 653 654 /* 655 Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 656 Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 657 */ 658 if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 659 /* construct the interpolation from the DMs */ 660 Mat p; 661 Vec rscale; 662 ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 663 dms[n-1] = pc->dm; 664 /* Separately create them so we do not get DMKSP interference between levels */ 665 for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 666 /* 667 Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level. 668 Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options. 669 But it is safe to use -dm_mat_type. 670 671 The mat type should not be hardcoded like this, we need to find a better way. 672 ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr); 673 */ 674 for (i=n-2; i>-1; i--) { 675 DMKSP kdm; 676 PetscBool dmhasrestrict, dmhasinject; 677 ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 678 if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 679 ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 680 /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 681 * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 682 kdm->ops->computerhs = NULL; 683 kdm->rhsctx = NULL; 684 if (!mglevels[i+1]->interpolate) { 685 ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 686 ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 687 if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 688 ierr = VecDestroy(&rscale);CHKERRQ(ierr); 689 ierr = MatDestroy(&p);CHKERRQ(ierr); 690 } 691 ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr); 692 if (dmhasrestrict && !mglevels[i+1]->restrct){ 693 ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr); 694 ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr); 695 ierr = MatDestroy(&p);CHKERRQ(ierr); 696 } 697 ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr); 698 if (dmhasinject && !mglevels[i+1]->inject){ 699 ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr); 700 ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr); 701 ierr = MatDestroy(&p);CHKERRQ(ierr); 702 } 703 } 704 705 for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 706 ierr = PetscFree(dms);CHKERRQ(ierr); 707 } 708 709 if (pc->dm && !pc->setupcalled) { 710 /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 711 ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 712 ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 713 } 714 715 if (mg->galerkin < PC_MG_GALERKIN_NONE) { 716 Mat A,B; 717 PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 718 MatReuse reuse = MAT_INITIAL_MATRIX; 719 720 if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 721 if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE; 722 if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 723 for (i=n-2; i>-1; i--) { 724 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"); 725 if (!mglevels[i+1]->interpolate) { 726 ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 727 } 728 if (!mglevels[i+1]->restrct) { 729 ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 730 } 731 if (reuse == MAT_REUSE_MATRIX) { 732 ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr); 733 } 734 if (doA) { 735 ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr); 736 } 737 if (doB) { 738 ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr); 739 } 740 /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 741 if (!doA && dAeqdB) { 742 if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);} 743 A = B; 744 } else if (!doA && reuse == MAT_INITIAL_MATRIX ) { 745 ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr); 746 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 747 } 748 if (!doB && dAeqdB) { 749 if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);} 750 B = A; 751 } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 752 ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 753 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 754 } 755 if (reuse == MAT_INITIAL_MATRIX) { 756 ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr); 757 ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr); 758 ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr); 759 } 760 dA = A; 761 dB = B; 762 } 763 } 764 if (needRestricts && pc->dm && pc->dm->x) { 765 /* need to restrict Jacobian location to coarser meshes for evaluation */ 766 for (i=n-2; i>-1; i--) { 767 Mat R; 768 Vec rscale; 769 if (!mglevels[i]->smoothd->dm->x) { 770 Vec *vecs; 771 ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 772 mglevels[i]->smoothd->dm->x = vecs[0]; 773 ierr = PetscFree(vecs);CHKERRQ(ierr); 774 } 775 ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 776 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 777 ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 778 ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 779 } 780 } 781 if (needRestricts && pc->dm) { 782 for (i=n-2; i>=0; i--) { 783 DM dmfine,dmcoarse; 784 Mat Restrict,Inject; 785 Vec rscale; 786 ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 787 ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 788 ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 789 ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 790 ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr); 791 ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 792 } 793 } 794 795 if (!pc->setupcalled) { 796 for (i=0; i<n; i++) { 797 ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 798 } 799 for (i=1; i<n; i++) { 800 if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 801 ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 802 } 803 } 804 /* insure that if either interpolation or restriction is set the other other one is set */ 805 for (i=1; i<n; i++) { 806 ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr); 807 ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr); 808 } 809 for (i=0; i<n-1; i++) { 810 if (!mglevels[i]->b) { 811 Vec *vec; 812 ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 813 ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 814 ierr = VecDestroy(vec);CHKERRQ(ierr); 815 ierr = PetscFree(vec);CHKERRQ(ierr); 816 } 817 if (!mglevels[i]->r && i) { 818 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 819 ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 820 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 821 } 822 if (!mglevels[i]->x) { 823 ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 824 ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 825 ierr = VecDestroy(&tvec);CHKERRQ(ierr); 826 } 827 } 828 if (n != 1 && !mglevels[n-1]->r) { 829 /* PCMGSetR() on the finest level if user did not supply it */ 830 Vec *vec; 831 ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 832 ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 833 ierr = VecDestroy(vec);CHKERRQ(ierr); 834 ierr = PetscFree(vec);CHKERRQ(ierr); 835 } 836 } 837 838 if (pc->dm) { 839 /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 840 for (i=0; i<n-1; i++) { 841 if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 842 } 843 } 844 845 for (i=1; i<n; i++) { 846 if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 847 /* if doing only down then initial guess is zero */ 848 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 849 } 850 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 851 ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 852 if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 853 pc->failedreason = PC_SUBPC_ERROR; 854 } 855 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 856 if (!mglevels[i]->residual) { 857 Mat mat; 858 ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 859 ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 860 } 861 } 862 for (i=1; i<n; i++) { 863 if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 864 Mat downmat,downpmat; 865 866 /* check if operators have been set for up, if not use down operators to set them */ 867 ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 868 if (!opsset) { 869 ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 870 ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 871 } 872 873 ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 874 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 875 ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 876 if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) { 877 pc->failedreason = PC_SUBPC_ERROR; 878 } 879 if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 880 } 881 } 882 883 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 884 ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 885 if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 886 pc->failedreason = PC_SUBPC_ERROR; 887 } 888 if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 889 890 /* 891 Dump the interpolation/restriction matrices plus the 892 Jacobian/stiffness on each level. This allows MATLAB users to 893 easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 894 895 Only support one or the other at the same time. 896 */ 897 #if defined(PETSC_USE_SOCKET_VIEWER) 898 ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 899 if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 900 dump = PETSC_FALSE; 901 #endif 902 ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 903 if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 904 905 if (viewer) { 906 for (i=1; i<n; i++) { 907 ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 908 } 909 for (i=0; i<n; i++) { 910 ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 911 ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 912 } 913 } 914 PetscFunctionReturn(0); 915 } 916 917 /* -------------------------------------------------------------------------------------*/ 918 919 PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 920 { 921 PC_MG *mg = (PC_MG *) pc->data; 922 923 PetscFunctionBegin; 924 *levels = mg->nlevels; 925 PetscFunctionReturn(0); 926 } 927 928 /*@ 929 PCMGGetLevels - Gets the number of levels to use with MG. 930 931 Not Collective 932 933 Input Parameter: 934 . pc - the preconditioner context 935 936 Output parameter: 937 . levels - the number of levels 938 939 Level: advanced 940 941 .keywords: MG, get, levels, multigrid 942 943 .seealso: PCMGSetLevels() 944 @*/ 945 PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 946 { 947 PetscErrorCode ierr; 948 949 PetscFunctionBegin; 950 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 951 PetscValidIntPointer(levels,2); 952 *levels = 0; 953 ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr); 954 PetscFunctionReturn(0); 955 } 956 957 /*@ 958 PCMGSetType - Determines the form of multigrid to use: 959 multiplicative, additive, full, or the Kaskade algorithm. 960 961 Logically Collective on PC 962 963 Input Parameters: 964 + pc - the preconditioner context 965 - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 966 PC_MG_FULL, PC_MG_KASKADE 967 968 Options Database Key: 969 . -pc_mg_type <form> - Sets <form>, one of multiplicative, 970 additive, full, kaskade 971 972 Level: advanced 973 974 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 975 976 .seealso: PCMGSetLevels() 977 @*/ 978 PetscErrorCode PCMGSetType(PC pc,PCMGType form) 979 { 980 PC_MG *mg = (PC_MG*)pc->data; 981 982 PetscFunctionBegin; 983 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 984 PetscValidLogicalCollectiveEnum(pc,form,2); 985 mg->am = form; 986 if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 987 else pc->ops->applyrichardson = NULL; 988 PetscFunctionReturn(0); 989 } 990 991 /*@ 992 PCMGGetType - Determines the form of multigrid to use: 993 multiplicative, additive, full, or the Kaskade algorithm. 994 995 Logically Collective on PC 996 997 Input Parameter: 998 . pc - the preconditioner context 999 1000 Output Parameter: 1001 . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 1002 1003 1004 Level: advanced 1005 1006 .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 1007 1008 .seealso: PCMGSetLevels() 1009 @*/ 1010 PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 1011 { 1012 PC_MG *mg = (PC_MG*)pc->data; 1013 1014 PetscFunctionBegin; 1015 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1016 *type = mg->am; 1017 PetscFunctionReturn(0); 1018 } 1019 1020 /*@ 1021 PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 1022 complicated cycling. 1023 1024 Logically Collective on PC 1025 1026 Input Parameters: 1027 + pc - the multigrid context 1028 - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 1029 1030 Options Database Key: 1031 . -pc_mg_cycle_type <v,w> - provide the cycle desired 1032 1033 Level: advanced 1034 1035 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 1036 1037 .seealso: PCMGSetCycleTypeOnLevel() 1038 @*/ 1039 PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 1040 { 1041 PC_MG *mg = (PC_MG*)pc->data; 1042 PC_MG_Levels **mglevels = mg->levels; 1043 PetscInt i,levels; 1044 1045 PetscFunctionBegin; 1046 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1047 PetscValidLogicalCollectiveEnum(pc,n,2); 1048 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1049 levels = mglevels[0]->levels; 1050 for (i=0; i<levels; i++) mglevels[i]->cycles = n; 1051 PetscFunctionReturn(0); 1052 } 1053 1054 /*@ 1055 PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 1056 of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 1057 1058 Logically Collective on PC 1059 1060 Input Parameters: 1061 + pc - the multigrid context 1062 - n - number of cycles (default is 1) 1063 1064 Options Database Key: 1065 . -pc_mg_multiplicative_cycles n 1066 1067 Level: advanced 1068 1069 Notes: 1070 This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 1071 1072 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 1073 1074 .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 1075 @*/ 1076 PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 1077 { 1078 PC_MG *mg = (PC_MG*)pc->data; 1079 1080 PetscFunctionBegin; 1081 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1082 PetscValidLogicalCollectiveInt(pc,n,2); 1083 mg->cyclesperpcapply = n; 1084 PetscFunctionReturn(0); 1085 } 1086 1087 PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1088 { 1089 PC_MG *mg = (PC_MG*)pc->data; 1090 1091 PetscFunctionBegin; 1092 mg->galerkin = use; 1093 PetscFunctionReturn(0); 1094 } 1095 1096 /*@ 1097 PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 1098 finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1099 1100 Logically Collective on PC 1101 1102 Input Parameters: 1103 + pc - the multigrid context 1104 - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1105 1106 Options Database Key: 1107 . -pc_mg_galerkin <both,pmat,mat,none> 1108 1109 Level: intermediate 1110 1111 Notes: 1112 Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 1113 use the PCMG construction of the coarser grids. 1114 1115 .keywords: MG, set, Galerkin 1116 1117 .seealso: PCMGGetGalerkin(), PCMGGalerkinType 1118 1119 @*/ 1120 PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1121 { 1122 PetscErrorCode ierr; 1123 1124 PetscFunctionBegin; 1125 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1126 ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr); 1127 PetscFunctionReturn(0); 1128 } 1129 1130 /*@ 1131 PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 1132 A_i-1 = r_i * A_i * p_i 1133 1134 Not Collective 1135 1136 Input Parameter: 1137 . pc - the multigrid context 1138 1139 Output Parameter: 1140 . galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL 1141 1142 Level: intermediate 1143 1144 .keywords: MG, set, Galerkin 1145 1146 .seealso: PCMGSetGalerkin(), PCMGGalerkinType 1147 1148 @*/ 1149 PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 1150 { 1151 PC_MG *mg = (PC_MG*)pc->data; 1152 1153 PetscFunctionBegin; 1154 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1155 *galerkin = mg->galerkin; 1156 PetscFunctionReturn(0); 1157 } 1158 1159 /*@ 1160 PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1161 on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1162 pre- and post-smoothing steps. 1163 1164 Logically Collective on PC 1165 1166 Input Parameters: 1167 + mg - the multigrid context 1168 - n - the number of smoothing steps 1169 1170 Options Database Key: 1171 + -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 1172 1173 Level: advanced 1174 1175 Notes: 1176 this does not set a value on the coarsest grid, since we assume that 1177 there is no separate smooth up on the coarsest grid. 1178 1179 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1180 1181 .seealso: PCMGSetDistinctSmoothUp() 1182 @*/ 1183 PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 1184 { 1185 PC_MG *mg = (PC_MG*)pc->data; 1186 PC_MG_Levels **mglevels = mg->levels; 1187 PetscErrorCode ierr; 1188 PetscInt i,levels; 1189 1190 PetscFunctionBegin; 1191 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1192 PetscValidLogicalCollectiveInt(pc,n,2); 1193 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1194 levels = mglevels[0]->levels; 1195 1196 for (i=1; i<levels; i++) { 1197 ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1198 ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 1199 mg->default_smoothu = n; 1200 mg->default_smoothd = n; 1201 } 1202 PetscFunctionReturn(0); 1203 } 1204 1205 /*@ 1206 PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels 1207 and adds the suffix _up to the options name 1208 1209 Logically Collective on PC 1210 1211 Input Parameters: 1212 . pc - the preconditioner context 1213 1214 Options Database Key: 1215 . -pc_mg_distinct_smoothup 1216 1217 Level: advanced 1218 1219 Notes: 1220 this does not set a value on the coarsest grid, since we assume that 1221 there is no separate smooth up on the coarsest grid. 1222 1223 .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1224 1225 .seealso: PCMGSetNumberSmooth() 1226 @*/ 1227 PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1228 { 1229 PC_MG *mg = (PC_MG*)pc->data; 1230 PC_MG_Levels **mglevels = mg->levels; 1231 PetscErrorCode ierr; 1232 PetscInt i,levels; 1233 KSP subksp; 1234 1235 PetscFunctionBegin; 1236 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1237 if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1238 levels = mglevels[0]->levels; 1239 1240 for (i=1; i<levels; i++) { 1241 const char *prefix = NULL; 1242 /* make sure smoother up and down are different */ 1243 ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr); 1244 ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr); 1245 ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr); 1246 ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr); 1247 } 1248 PetscFunctionReturn(0); 1249 } 1250 1251 /* ----------------------------------------------------------------------------------------*/ 1252 1253 /*MC 1254 PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 1255 information about the coarser grid matrices and restriction/interpolation operators. 1256 1257 Options Database Keys: 1258 + -pc_mg_levels <nlevels> - number of levels including finest 1259 . -pc_mg_cycle_type <v,w> - provide the cycle desired 1260 . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 1261 . -pc_mg_log - log information about time spent on each level of the solver 1262 . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 1263 . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 1264 . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 1265 . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1266 to the Socket viewer for reading from MATLAB. 1267 - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 1268 to the binary output file called binaryoutput 1269 1270 Notes: 1271 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 1272 1273 When run with a single level the smoother options are used on that level NOT the coarse grid solver options 1274 1275 When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 1276 is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 1277 (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 1278 residual is computed at the end of each cycle. 1279 1280 Level: intermediate 1281 1282 Concepts: multigrid/multilevel 1283 1284 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1285 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), 1286 PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 1287 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 1288 PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 1289 M*/ 1290 1291 PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 1292 { 1293 PC_MG *mg; 1294 PetscErrorCode ierr; 1295 1296 PetscFunctionBegin; 1297 ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1298 pc->data = (void*)mg; 1299 mg->nlevels = -1; 1300 mg->am = PC_MG_MULTIPLICATIVE; 1301 mg->galerkin = PC_MG_GALERKIN_NONE; 1302 1303 pc->useAmat = PETSC_TRUE; 1304 1305 pc->ops->apply = PCApply_MG; 1306 pc->ops->setup = PCSetUp_MG; 1307 pc->ops->reset = PCReset_MG; 1308 pc->ops->destroy = PCDestroy_MG; 1309 pc->ops->setfromoptions = PCSetFromOptions_MG; 1310 pc->ops->view = PCView_MG; 1311 1312 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 1313 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr); 1314 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr); 1315 PetscFunctionReturn(0); 1316 } 1317