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