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