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