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