1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 34b9ad928SBarry Smith /* 44b9ad928SBarry Smith Defines the multigrid preconditioner interface. 54b9ad928SBarry Smith */ 67c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 74b9ad928SBarry Smith 84b9ad928SBarry Smith 94b9ad928SBarry Smith #undef __FUNCT__ 109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private" 1131567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 124b9ad928SBarry Smith { 1331567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1431567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 156849ba73SBarry Smith PetscErrorCode ierr; 1631567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 174b9ad928SBarry Smith 184b9ad928SBarry Smith PetscFunctionBegin; 194b9ad928SBarry Smith 2032cf1786SBarry Smith if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2131567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 2232cf1786SBarry Smith if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2331567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 2432cf1786SBarry Smith if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);} 2531567311SBarry Smith ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 2632cf1786SBarry Smith if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);} 274b9ad928SBarry Smith 284b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 2931567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 304b9ad928SBarry Smith PetscReal rnorm; 3131567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 324b9ad928SBarry Smith if (rnorm <= mg->ttol) { 3370441072SBarry Smith if (rnorm < mg->abstol) { 344d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 351e2582c4SBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr); 364b9ad928SBarry Smith } else { 374d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 381e2582c4SBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than relative tolerance times initial residual norm %G\n",rnorm,mg->ttol);CHKERRQ(ierr); 394b9ad928SBarry Smith } 404b9ad928SBarry Smith PetscFunctionReturn(0); 414b9ad928SBarry Smith } 424b9ad928SBarry Smith } 434b9ad928SBarry Smith 4431567311SBarry Smith mgc = *(mglevelsin - 1); 4532cf1786SBarry Smith if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 4631567311SBarry Smith ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 4732cf1786SBarry Smith if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 48efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 494b9ad928SBarry Smith while (cycles--) { 5031567311SBarry Smith ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 514b9ad928SBarry Smith } 5232cf1786SBarry Smith if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5331567311SBarry Smith ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 5432cf1786SBarry Smith if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5532cf1786SBarry Smith if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 5631567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 5732cf1786SBarry Smith if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 584b9ad928SBarry Smith } 594b9ad928SBarry Smith PetscFunctionReturn(0); 604b9ad928SBarry Smith } 614b9ad928SBarry Smith 624b9ad928SBarry Smith #undef __FUNCT__ 634b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG" 647319c654SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscTruth zeroguess,PetscInt *outits,PCRichardsonConvergedReason *reason) 654b9ad928SBarry Smith { 66f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 67f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 68dfbe8321SBarry Smith PetscErrorCode ierr; 69f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 704b9ad928SBarry Smith 714b9ad928SBarry Smith PetscFunctionBegin; 72f3fbd535SBarry Smith mglevels[levels-1]->b = b; 73f3fbd535SBarry Smith mglevels[levels-1]->x = x; 744b9ad928SBarry Smith 7531567311SBarry Smith mg->rtol = rtol; 7631567311SBarry Smith mg->abstol = abstol; 7731567311SBarry Smith mg->dtol = dtol; 784b9ad928SBarry Smith if (rtol) { 794b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 804b9ad928SBarry Smith PetscReal rnorm; 817319c654SBarry Smith if (zeroguess) { 827319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 837319c654SBarry Smith } else { 84f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 854b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 867319c654SBarry Smith } 8731567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 8870441072SBarry Smith } else if (abstol) { 8931567311SBarry Smith mg->ttol = abstol; 904b9ad928SBarry Smith } else { 9131567311SBarry Smith mg->ttol = 0.0; 924b9ad928SBarry Smith } 934b9ad928SBarry Smith 944d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 954d0a8057SBarry Smith stop prematurely do to small residual */ 964d0a8057SBarry Smith for (i=1; i<levels; i++) { 97f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 98f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 99f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 1004b9ad928SBarry Smith } 1014d0a8057SBarry Smith } 1024d0a8057SBarry Smith 1034d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1044d0a8057SBarry Smith for (i=0; i<its; i++) { 105f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1064d0a8057SBarry Smith if (*reason) break; 1074d0a8057SBarry Smith } 1084d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1094d0a8057SBarry Smith *outits = i; 1104b9ad928SBarry Smith PetscFunctionReturn(0); 1114b9ad928SBarry Smith } 1124b9ad928SBarry Smith 1134b9ad928SBarry Smith #undef __FUNCT__ 1149dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 1154b9ad928SBarry Smith /*@C 11697177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 1174b9ad928SBarry Smith Must be called before any other MG routine. 1184b9ad928SBarry Smith 1194b9ad928SBarry Smith Collective on PC 1204b9ad928SBarry Smith 1214b9ad928SBarry Smith Input Parameters: 1224b9ad928SBarry Smith + pc - the preconditioner context 1234b9ad928SBarry Smith . levels - the number of levels 1244b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 1254b9ad928SBarry Smith on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 1264b9ad928SBarry Smith 1274b9ad928SBarry Smith Level: intermediate 1284b9ad928SBarry Smith 1294b9ad928SBarry Smith Notes: 1304b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 1314b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 1324b9ad928SBarry Smith 1334b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 1344b9ad928SBarry Smith 13597177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 1364b9ad928SBarry Smith @*/ 13797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 1384b9ad928SBarry Smith { 139dfbe8321SBarry Smith PetscErrorCode ierr; 140f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 141f3fbd535SBarry Smith MPI_Comm comm = ((PetscObject)pc)->comm; 142f3fbd535SBarry Smith PC_MG_Levels **mglevels; 143f3fbd535SBarry Smith PetscInt i; 144f3fbd535SBarry Smith PetscMPIInt size; 145f3fbd535SBarry Smith const char *prefix; 146f3fbd535SBarry Smith PC ipc; 1474b9ad928SBarry Smith 1484b9ad928SBarry Smith PetscFunctionBegin; 1494482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 150f3fbd535SBarry Smith if (mg->nlevels > -1) { 151f3fbd535SBarry Smith SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n make sure that you call PCMGSetLevels() before KSPSetFromOptions()"); 1524b9ad928SBarry Smith } 153f3fbd535SBarry Smith if (mg->levels) SETERRQ(PETSC_ERR_PLIB,"Internal error in PETSc, this array should not yet exist"); 154f3fbd535SBarry Smith 155f3fbd535SBarry Smith mg->nlevels = levels; 156f3fbd535SBarry Smith 157f3fbd535SBarry Smith ierr = PetscMalloc(levels*sizeof(PC_MG*),&mglevels);CHKERRQ(ierr); 158f3fbd535SBarry Smith ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 159f3fbd535SBarry Smith 160f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 161f3fbd535SBarry Smith 162f3fbd535SBarry Smith for (i=0; i<levels; i++) { 163f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG_Levels,&mglevels[i]);CHKERRQ(ierr); 164f3fbd535SBarry Smith mglevels[i]->level = i; 165f3fbd535SBarry Smith mglevels[i]->levels = levels; 166f3fbd535SBarry Smith mglevels[i]->cycles = PC_MG_CYCLE_V; 167f3fbd535SBarry Smith mglevels[i]->galerkin = PETSC_FALSE; 168f3fbd535SBarry Smith mglevels[i]->galerkinused = PETSC_FALSE; 16931567311SBarry Smith mg->default_smoothu = 1; 17031567311SBarry Smith mg->default_smoothd = 1; 171f3fbd535SBarry Smith 172f3fbd535SBarry Smith if (comms) comm = comms[i]; 173f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 174f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 17531567311SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 176f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 177f3fbd535SBarry Smith 178f3fbd535SBarry Smith /* do special stuff for coarse grid */ 179f3fbd535SBarry Smith if (!i && levels > 1) { 180f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 181f3fbd535SBarry Smith 182f3fbd535SBarry Smith /* coarse solve is (redundant) LU by default */ 183f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 184f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 185f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 186f3fbd535SBarry Smith if (size > 1) { 187f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 188f3fbd535SBarry Smith } else { 189f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 190f3fbd535SBarry Smith } 191f3fbd535SBarry Smith 192f3fbd535SBarry Smith } else { 193f3fbd535SBarry Smith char tprefix[128]; 194f3fbd535SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 195f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 196f3fbd535SBarry Smith } 197f3fbd535SBarry Smith ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr); 198f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 19931567311SBarry Smith mg->rtol = 0.0; 20031567311SBarry Smith mg->abstol = 0.0; 20131567311SBarry Smith mg->dtol = 0.0; 20231567311SBarry Smith mg->ttol = 0.0; 20331567311SBarry Smith mg->eventsmoothsetup = 0; 20431567311SBarry Smith mg->eventsmoothsolve = 0; 20531567311SBarry Smith mg->eventresidual = 0; 20631567311SBarry Smith mg->eventinterprestrict = 0; 20731567311SBarry Smith mg->cyclesperpcapply = 1; 208f3fbd535SBarry Smith } 20931567311SBarry Smith mg->am = PC_MG_MULTIPLICATIVE; 210f3fbd535SBarry Smith mg->levels = mglevels; 2114b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_MG; 2124b9ad928SBarry Smith PetscFunctionReturn(0); 2134b9ad928SBarry Smith } 2144b9ad928SBarry Smith 2154b9ad928SBarry Smith #undef __FUNCT__ 216*c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG_Private" 217*c07bf074SBarry Smith PetscErrorCode PCDestroy_MG_Private(PC pc) 218f3fbd535SBarry Smith { 219f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 220f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 221f3fbd535SBarry Smith PetscErrorCode ierr; 222f3fbd535SBarry Smith PetscInt i,n; 223f3fbd535SBarry Smith 224f3fbd535SBarry Smith PetscFunctionBegin; 225f3fbd535SBarry Smith if (mglevels) { 226f3fbd535SBarry Smith n = mglevels[0]->levels; 227f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 228f3fbd535SBarry Smith if (mglevels[i+1]->r) {ierr = VecDestroy(mglevels[i+1]->r);CHKERRQ(ierr);} 229f3fbd535SBarry Smith if (mglevels[i]->b) {ierr = VecDestroy(mglevels[i]->b);CHKERRQ(ierr);} 230f3fbd535SBarry Smith if (mglevels[i]->x) {ierr = VecDestroy(mglevels[i]->x);CHKERRQ(ierr);} 231f3fbd535SBarry Smith if (mglevels[i+1]->restrct) {ierr = MatDestroy(mglevels[i+1]->restrct);CHKERRQ(ierr);} 232f3fbd535SBarry Smith if (mglevels[i+1]->interpolate) {ierr = MatDestroy(mglevels[i+1]->interpolate);CHKERRQ(ierr);} 233f3fbd535SBarry Smith } 234f3fbd535SBarry Smith 235f3fbd535SBarry Smith for (i=0; i<n; i++) { 236f3fbd535SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 237f3fbd535SBarry Smith ierr = KSPDestroy(mglevels[i]->smoothd);CHKERRQ(ierr); 238f3fbd535SBarry Smith } 239f3fbd535SBarry Smith ierr = KSPDestroy(mglevels[i]->smoothu);CHKERRQ(ierr); 240f3fbd535SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 241f3fbd535SBarry Smith } 242f3fbd535SBarry Smith ierr = PetscFree(mglevels);CHKERRQ(ierr); 243f3fbd535SBarry Smith } 244*c07bf074SBarry Smith mg->nlevels = -1; 245*c07bf074SBarry Smith mg->levels = PETSC_NULL; 246*c07bf074SBarry Smith PetscFunctionReturn(0); 247*c07bf074SBarry Smith } 248*c07bf074SBarry Smith 249*c07bf074SBarry Smith #undef __FUNCT__ 250*c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG" 251*c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 252*c07bf074SBarry Smith { 253*c07bf074SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 254*c07bf074SBarry Smith PetscErrorCode ierr; 255*c07bf074SBarry Smith 256*c07bf074SBarry Smith PetscFunctionBegin; 257*c07bf074SBarry Smith ierr = PCDestroy_MG_Private(pc);CHKERRQ(ierr); 258f3fbd535SBarry Smith ierr = PetscFree(mg);CHKERRQ(ierr); 259f3fbd535SBarry Smith PetscFunctionReturn(0); 260f3fbd535SBarry Smith } 261f3fbd535SBarry Smith 262f3fbd535SBarry Smith 263f3fbd535SBarry Smith 26431567311SBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 265f3fbd535SBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 26631567311SBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 267f3fbd535SBarry Smith 268f3fbd535SBarry Smith /* 269f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 270f3fbd535SBarry Smith or full cycle of multigrid. 271f3fbd535SBarry Smith 272f3fbd535SBarry Smith Note: 273f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 274f3fbd535SBarry Smith */ 275f3fbd535SBarry Smith #undef __FUNCT__ 276f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG" 277f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 278f3fbd535SBarry Smith { 279f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 280f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 281f3fbd535SBarry Smith PetscErrorCode ierr; 282f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 283f3fbd535SBarry Smith 284f3fbd535SBarry Smith PetscFunctionBegin; 285f3fbd535SBarry Smith mglevels[levels-1]->b = b; 286f3fbd535SBarry Smith mglevels[levels-1]->x = x; 28731567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 288f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 28931567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 290f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr); 291f3fbd535SBarry Smith } 292f3fbd535SBarry Smith } 29331567311SBarry Smith else if (mg->am == PC_MG_ADDITIVE) { 29431567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 295f3fbd535SBarry Smith } 29631567311SBarry Smith else if (mg->am == PC_MG_KASKADE) { 29731567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 298f3fbd535SBarry Smith } 299f3fbd535SBarry Smith else { 300f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 301f3fbd535SBarry Smith } 302f3fbd535SBarry Smith PetscFunctionReturn(0); 303f3fbd535SBarry Smith } 304f3fbd535SBarry Smith 305f3fbd535SBarry Smith 306f3fbd535SBarry Smith #undef __FUNCT__ 307f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 308f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc) 309f3fbd535SBarry Smith { 310f3fbd535SBarry Smith PetscErrorCode ierr; 311f3fbd535SBarry Smith PetscInt m,levels = 1,cycles; 312f3fbd535SBarry Smith PetscTruth flg; 313f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 314f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 315f3fbd535SBarry Smith PCMGType mgtype; 316f3fbd535SBarry Smith PCMGCycleType mgctype; 317f3fbd535SBarry Smith 318f3fbd535SBarry Smith PetscFunctionBegin; 319f3fbd535SBarry Smith ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 320f3fbd535SBarry Smith if (!pc->data) { 321f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 322f3fbd535SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 323f3fbd535SBarry Smith mglevels = mg->levels; 324f3fbd535SBarry Smith } 325f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 326f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 327f3fbd535SBarry Smith if (flg) { 328f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 329f3fbd535SBarry Smith }; 330f3fbd535SBarry Smith flg = PETSC_FALSE; 331f3fbd535SBarry Smith ierr = PetscOptionsTruth("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 332f3fbd535SBarry Smith if (flg) { 333f3fbd535SBarry Smith ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 334f3fbd535SBarry Smith } 335f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 336f3fbd535SBarry Smith if (flg) { 337f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 338f3fbd535SBarry Smith } 339f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 340f3fbd535SBarry Smith if (flg) { 341f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 342f3fbd535SBarry Smith } 34331567311SBarry Smith mgtype = mg->am; 344f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 345f3fbd535SBarry Smith if (flg) { 346f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 347f3fbd535SBarry Smith } 34831567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 34931567311SBarry Smith ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 350f3fbd535SBarry Smith if (flg) { 351f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 352f3fbd535SBarry Smith } 353f3fbd535SBarry Smith } 354f3fbd535SBarry Smith flg = PETSC_FALSE; 355f3fbd535SBarry Smith ierr = PetscOptionsTruth("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 356f3fbd535SBarry Smith if (flg) { 357f3fbd535SBarry Smith PetscInt i; 358f3fbd535SBarry Smith char eventname[128]; 359f3fbd535SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 360f3fbd535SBarry Smith levels = mglevels[0]->levels; 361f3fbd535SBarry Smith for (i=0; i<levels; i++) { 362f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 36331567311SBarry Smith ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg->eventsmoothsetup);CHKERRQ(ierr); 364f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 36531567311SBarry Smith ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg->eventsmoothsolve);CHKERRQ(ierr); 366f3fbd535SBarry Smith if (i) { 367f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 36831567311SBarry Smith ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg->eventresidual);CHKERRQ(ierr); 369f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 37031567311SBarry Smith ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg->eventinterprestrict);CHKERRQ(ierr); 371f3fbd535SBarry Smith } 372f3fbd535SBarry Smith } 373f3fbd535SBarry Smith } 374f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 375f3fbd535SBarry Smith PetscFunctionReturn(0); 376f3fbd535SBarry Smith } 377f3fbd535SBarry Smith 378f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 379f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 380f3fbd535SBarry Smith 381f3fbd535SBarry Smith #undef __FUNCT__ 382f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG" 383f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 384f3fbd535SBarry Smith { 385f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 386f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 387f3fbd535SBarry Smith PetscErrorCode ierr; 388f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 389f3fbd535SBarry Smith PetscTruth iascii; 390f3fbd535SBarry Smith 391f3fbd535SBarry Smith PetscFunctionBegin; 392f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 393f3fbd535SBarry Smith if (iascii) { 39431567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,(mglevels[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr); 39531567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 39631567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 397f3fbd535SBarry Smith } 398f3fbd535SBarry Smith if (mglevels[0]->galerkin) { 399f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 400f3fbd535SBarry Smith } 401f3fbd535SBarry Smith for (i=0; i<levels; i++) { 402f3fbd535SBarry Smith if (!i) { 40331567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D presmooths=%D postsmooths=%D -----\n",i,mg->default_smoothd,mg->default_smoothu);CHKERRQ(ierr); 404f3fbd535SBarry Smith } else { 40531567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr); 406f3fbd535SBarry Smith } 407f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 408f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 409f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 410f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 411f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 412f3fbd535SBarry Smith } else if (i){ 41331567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr); 414f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 415f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 416f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 417f3fbd535SBarry Smith } 418f3fbd535SBarry Smith } 419f3fbd535SBarry Smith } else { 420f3fbd535SBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 421f3fbd535SBarry Smith } 422f3fbd535SBarry Smith PetscFunctionReturn(0); 423f3fbd535SBarry Smith } 424f3fbd535SBarry Smith 425f3fbd535SBarry Smith /* 426f3fbd535SBarry Smith Calls setup for the KSP on each level 427f3fbd535SBarry Smith */ 428f3fbd535SBarry Smith #undef __FUNCT__ 429f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG" 430f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 431f3fbd535SBarry Smith { 432f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 433f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 434f3fbd535SBarry Smith PetscErrorCode ierr; 435f3fbd535SBarry Smith PetscInt i,n = mglevels[0]->levels; 436f3fbd535SBarry Smith PC cpc,mpc; 437f3fbd535SBarry Smith PetscTruth preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset; 438f3fbd535SBarry Smith PetscViewerASCIIMonitor ascii; 439f3fbd535SBarry Smith PetscViewer viewer = PETSC_NULL; 440f3fbd535SBarry Smith MPI_Comm comm; 441f3fbd535SBarry Smith Mat dA,dB; 442f3fbd535SBarry Smith MatStructure uflag; 443f3fbd535SBarry Smith Vec tvec; 444f3fbd535SBarry Smith 445f3fbd535SBarry Smith PetscFunctionBegin; 446f3fbd535SBarry Smith 447f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 448f3fbd535SBarry Smith /* so use those from global PC */ 449f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 450f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 451f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 452f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr); 4531338a6b9SJed Brown if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) { 454f3fbd535SBarry Smith 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); 455f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 456f3fbd535SBarry Smith } 457f3fbd535SBarry Smith 458f3fbd535SBarry Smith if (mglevels[0]->galerkin) { 459f3fbd535SBarry Smith Mat B; 460f3fbd535SBarry Smith mglevels[0]->galerkinused = PETSC_TRUE; 461f3fbd535SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 462f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 463f3fbd535SBarry Smith if (!pc->setupcalled) { 464f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 465f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 466f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 467f3fbd535SBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 468f3fbd535SBarry Smith dB = B; 469f3fbd535SBarry Smith } 470f3fbd535SBarry Smith ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr); 471f3fbd535SBarry Smith } else { 472f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 473f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 474f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 475f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 476f3fbd535SBarry Smith dB = B; 477f3fbd535SBarry Smith } 478f3fbd535SBarry Smith } 479f3fbd535SBarry Smith } 480f3fbd535SBarry Smith 481f3fbd535SBarry Smith if (!pc->setupcalled) { 482f3fbd535SBarry Smith ierr = PetscOptionsGetTruth(0,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr); 483f3fbd535SBarry Smith 484f3fbd535SBarry Smith for (i=0; i<n; i++) { 485f3fbd535SBarry Smith if (monitor) { 486f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr); 487f3fbd535SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 488f3fbd535SBarry Smith ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 489f3fbd535SBarry Smith } 490f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 491f3fbd535SBarry Smith } 492f3fbd535SBarry Smith for (i=1; i<n; i++) { 493f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 494f3fbd535SBarry Smith if (monitor) { 495f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr); 496f3fbd535SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 497f3fbd535SBarry Smith ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 498f3fbd535SBarry Smith } 499f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 500f3fbd535SBarry Smith } 501f3fbd535SBarry Smith } 502f3fbd535SBarry Smith for (i=1; i<n; i++) { 503f3fbd535SBarry Smith if (!mglevels[i]->residual) { 504f3fbd535SBarry Smith Mat mat; 505f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 506f3fbd535SBarry Smith ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 507f3fbd535SBarry Smith } 508f3fbd535SBarry Smith if (mglevels[i]->restrct && !mglevels[i]->interpolate) { 509f3fbd535SBarry Smith ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr); 510f3fbd535SBarry Smith } 511f3fbd535SBarry Smith if (!mglevels[i]->restrct && mglevels[i]->interpolate) { 512f3fbd535SBarry Smith ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr); 513f3fbd535SBarry Smith } 514f3fbd535SBarry Smith #if defined(PETSC_USE_DEBUG) 515f3fbd535SBarry Smith if (!mglevels[i]->restrct || !mglevels[i]->interpolate) { 516f3fbd535SBarry Smith SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 517f3fbd535SBarry Smith } 518f3fbd535SBarry Smith #endif 519f3fbd535SBarry Smith } 520f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 521f3fbd535SBarry Smith if (!mglevels[i]->b) { 522f3fbd535SBarry Smith Vec *vec; 523f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 524f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 525f3fbd535SBarry Smith ierr = VecDestroy(*vec);CHKERRQ(ierr); 526f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 527f3fbd535SBarry Smith } 528f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 529f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 530f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 531f3fbd535SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 532f3fbd535SBarry Smith } 533f3fbd535SBarry Smith if (!mglevels[i]->x) { 534f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 535f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 536f3fbd535SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 537f3fbd535SBarry Smith } 538f3fbd535SBarry Smith } 539f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 540f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 541f3fbd535SBarry Smith Vec *vec; 542f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 543f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 544f3fbd535SBarry Smith ierr = VecDestroy(*vec);CHKERRQ(ierr); 545f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 546f3fbd535SBarry Smith } 547f3fbd535SBarry Smith } 548f3fbd535SBarry Smith 549f3fbd535SBarry Smith 550f3fbd535SBarry Smith for (i=1; i<n; i++) { 551f3fbd535SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 552f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 553f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 554f3fbd535SBarry Smith } 55531567311SBarry Smith if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 556f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 55731567311SBarry Smith if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 558f3fbd535SBarry Smith } 559f3fbd535SBarry Smith for (i=1; i<n; i++) { 560f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 561f3fbd535SBarry Smith Mat downmat,downpmat; 562f3fbd535SBarry Smith MatStructure matflag; 563f3fbd535SBarry Smith PetscTruth opsset; 564f3fbd535SBarry Smith 565f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 566f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 567f3fbd535SBarry Smith if (!opsset) { 568f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 569f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 570f3fbd535SBarry Smith } 571f3fbd535SBarry Smith 572f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 57331567311SBarry Smith if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 574f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 57531567311SBarry Smith if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 576f3fbd535SBarry Smith } 577f3fbd535SBarry Smith } 578f3fbd535SBarry Smith 579f3fbd535SBarry Smith /* 580f3fbd535SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 581f3fbd535SBarry Smith */ 582f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 583f3fbd535SBarry Smith if (preonly) { 584f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 585f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 586f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 587f3fbd535SBarry Smith if (!lu && !redundant && !cholesky) { 588f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 589f3fbd535SBarry Smith } 590f3fbd535SBarry Smith } 591f3fbd535SBarry Smith 592f3fbd535SBarry Smith if (!pc->setupcalled) { 593f3fbd535SBarry Smith if (monitor) { 594f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr); 595f3fbd535SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 596f3fbd535SBarry Smith ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 597f3fbd535SBarry Smith } 598f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 599f3fbd535SBarry Smith } 600f3fbd535SBarry Smith 60131567311SBarry Smith if (mg->eventsmoothsetup) {ierr = PetscLogEventBegin(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 602f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 60331567311SBarry Smith if (mg->eventsmoothsetup) {ierr = PetscLogEventEnd(mg->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 604f3fbd535SBarry Smith 605f3fbd535SBarry Smith /* 606f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 607f3fbd535SBarry Smith Jacobian/stiffness on each level. This allows Matlab users to 608f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 609f3fbd535SBarry Smith 610f3fbd535SBarry Smith Only support one or the other at the same time. 611f3fbd535SBarry Smith */ 612f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 613f3fbd535SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 614f3fbd535SBarry Smith if (dump) { 615f3fbd535SBarry Smith viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 616f3fbd535SBarry Smith } 617f3fbd535SBarry Smith dump = PETSC_FALSE; 618f3fbd535SBarry Smith #endif 619f3fbd535SBarry Smith ierr = PetscOptionsGetTruth(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 620f3fbd535SBarry Smith if (dump) { 621f3fbd535SBarry Smith viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 622f3fbd535SBarry Smith } 623f3fbd535SBarry Smith 624f3fbd535SBarry Smith if (viewer) { 625f3fbd535SBarry Smith for (i=1; i<n; i++) { 626f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 627f3fbd535SBarry Smith } 628f3fbd535SBarry Smith for (i=0; i<n; i++) { 629f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 630f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 631f3fbd535SBarry Smith } 632f3fbd535SBarry Smith } 633f3fbd535SBarry Smith PetscFunctionReturn(0); 634f3fbd535SBarry Smith } 635f3fbd535SBarry Smith 636f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 637f3fbd535SBarry Smith 638f3fbd535SBarry Smith #undef __FUNCT__ 6399dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 6404b9ad928SBarry Smith /*@ 64197177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 6424b9ad928SBarry Smith 6434b9ad928SBarry Smith Not Collective 6444b9ad928SBarry Smith 6454b9ad928SBarry Smith Input Parameter: 6464b9ad928SBarry Smith . pc - the preconditioner context 6474b9ad928SBarry Smith 6484b9ad928SBarry Smith Output parameter: 6494b9ad928SBarry Smith . levels - the number of levels 6504b9ad928SBarry Smith 6514b9ad928SBarry Smith Level: advanced 6524b9ad928SBarry Smith 6534b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 6544b9ad928SBarry Smith 65597177400SBarry Smith .seealso: PCMGSetLevels() 6564b9ad928SBarry Smith @*/ 65797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels) 6584b9ad928SBarry Smith { 659f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 6604b9ad928SBarry Smith 6614b9ad928SBarry Smith PetscFunctionBegin; 6624482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6634482741eSBarry Smith PetscValidIntPointer(levels,2); 664f3fbd535SBarry Smith *levels = mg->nlevels; 6654b9ad928SBarry Smith PetscFunctionReturn(0); 6664b9ad928SBarry Smith } 6674b9ad928SBarry Smith 6684b9ad928SBarry Smith #undef __FUNCT__ 6699dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 6704b9ad928SBarry Smith /*@ 67197177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 6724b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 6734b9ad928SBarry Smith 6744b9ad928SBarry Smith Collective on PC 6754b9ad928SBarry Smith 6764b9ad928SBarry Smith Input Parameters: 6774b9ad928SBarry Smith + pc - the preconditioner context 6789dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 6799dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 6804b9ad928SBarry Smith 6814b9ad928SBarry Smith Options Database Key: 6824b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 6834b9ad928SBarry Smith additive, full, kaskade 6844b9ad928SBarry Smith 6854b9ad928SBarry Smith Level: advanced 6864b9ad928SBarry Smith 6874b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 6884b9ad928SBarry Smith 68997177400SBarry Smith .seealso: PCMGSetLevels() 6904b9ad928SBarry Smith @*/ 6919dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form) 6924b9ad928SBarry Smith { 693f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 6944b9ad928SBarry Smith 6954b9ad928SBarry Smith PetscFunctionBegin; 6964482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 69731567311SBarry Smith mg->am = form; 6989dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 6994b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 7004b9ad928SBarry Smith PetscFunctionReturn(0); 7014b9ad928SBarry Smith } 7024b9ad928SBarry Smith 7034b9ad928SBarry Smith #undef __FUNCT__ 7040d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 7054b9ad928SBarry Smith /*@ 7060d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 7074b9ad928SBarry Smith complicated cycling. 7084b9ad928SBarry Smith 7094b9ad928SBarry Smith Collective on PC 7104b9ad928SBarry Smith 7114b9ad928SBarry Smith Input Parameters: 712c2be2410SBarry Smith + pc - the multigrid context 7130d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 7144b9ad928SBarry Smith 7154b9ad928SBarry Smith Options Database Key: 7160d353602SBarry Smith $ -pc_mg_cycle_type v or w 7174b9ad928SBarry Smith 7184b9ad928SBarry Smith Level: advanced 7194b9ad928SBarry Smith 7204b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 7214b9ad928SBarry Smith 7220d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 7234b9ad928SBarry Smith @*/ 7240d353602SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n) 7254b9ad928SBarry Smith { 726f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 727f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 72879416396SBarry Smith PetscInt i,levels; 7294b9ad928SBarry Smith 7304b9ad928SBarry Smith PetscFunctionBegin; 7314482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 732f3fbd535SBarry Smith if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 733f3fbd535SBarry Smith levels = mglevels[0]->levels; 7344b9ad928SBarry Smith 7354b9ad928SBarry Smith for (i=0; i<levels; i++) { 736f3fbd535SBarry Smith mglevels[i]->cycles = n; 7374b9ad928SBarry Smith } 7384b9ad928SBarry Smith PetscFunctionReturn(0); 7394b9ad928SBarry Smith } 7404b9ad928SBarry Smith 7414b9ad928SBarry Smith #undef __FUNCT__ 7428cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 7438cc2d5dfSBarry Smith /*@ 7448cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 7458cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 7468cc2d5dfSBarry Smith 7478cc2d5dfSBarry Smith Collective on PC 7488cc2d5dfSBarry Smith 7498cc2d5dfSBarry Smith Input Parameters: 7508cc2d5dfSBarry Smith + pc - the multigrid context 7518cc2d5dfSBarry Smith - n - number of cycles (default is 1) 7528cc2d5dfSBarry Smith 7538cc2d5dfSBarry Smith Options Database Key: 7548cc2d5dfSBarry Smith $ -pc_mg_multiplicative_cycles n 7558cc2d5dfSBarry Smith 7568cc2d5dfSBarry Smith Level: advanced 7578cc2d5dfSBarry Smith 7588cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 7598cc2d5dfSBarry Smith 7608cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 7618cc2d5dfSBarry Smith 7628cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 7638cc2d5dfSBarry Smith @*/ 7648cc2d5dfSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 7658cc2d5dfSBarry Smith { 766f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 767f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 7688cc2d5dfSBarry Smith PetscInt i,levels; 7698cc2d5dfSBarry Smith 7708cc2d5dfSBarry Smith PetscFunctionBegin; 7718cc2d5dfSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 772f3fbd535SBarry Smith if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 773f3fbd535SBarry Smith levels = mglevels[0]->levels; 7748cc2d5dfSBarry Smith 7758cc2d5dfSBarry Smith for (i=0; i<levels; i++) { 77631567311SBarry Smith mg->cyclesperpcapply = n; 7778cc2d5dfSBarry Smith } 7788cc2d5dfSBarry Smith PetscFunctionReturn(0); 7798cc2d5dfSBarry Smith } 7808cc2d5dfSBarry Smith 7818cc2d5dfSBarry Smith #undef __FUNCT__ 7829dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 783c2be2410SBarry Smith /*@ 78497177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 785c2be2410SBarry Smith finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 786c2be2410SBarry Smith 787c2be2410SBarry Smith Collective on PC 788c2be2410SBarry Smith 789c2be2410SBarry Smith Input Parameters: 7903fc8bf9cSBarry Smith . pc - the multigrid context 791c2be2410SBarry Smith 792c2be2410SBarry Smith Options Database Key: 793c2be2410SBarry Smith $ -pc_mg_galerkin 794c2be2410SBarry Smith 795c2be2410SBarry Smith Level: intermediate 796c2be2410SBarry Smith 797c2be2410SBarry Smith .keywords: MG, set, Galerkin 798c2be2410SBarry Smith 7993fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 8003fc8bf9cSBarry Smith 801c2be2410SBarry Smith @*/ 80297177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 803c2be2410SBarry Smith { 804f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 805f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 806c2be2410SBarry Smith PetscInt i,levels; 807c2be2410SBarry Smith 808c2be2410SBarry Smith PetscFunctionBegin; 809c2be2410SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 810f3fbd535SBarry Smith if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 811f3fbd535SBarry Smith levels = mglevels[0]->levels; 812c2be2410SBarry Smith 813c2be2410SBarry Smith for (i=0; i<levels; i++) { 814f3fbd535SBarry Smith mglevels[i]->galerkin = PETSC_TRUE; 815c2be2410SBarry Smith } 816c2be2410SBarry Smith PetscFunctionReturn(0); 817c2be2410SBarry Smith } 818c2be2410SBarry Smith 819c2be2410SBarry Smith #undef __FUNCT__ 8203fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 8213fc8bf9cSBarry Smith /*@ 8223fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 8233fc8bf9cSBarry Smith A_i-1 = r_i * A_i * r_i^t 8243fc8bf9cSBarry Smith 8253fc8bf9cSBarry Smith Not Collective 8263fc8bf9cSBarry Smith 8273fc8bf9cSBarry Smith Input Parameter: 8283fc8bf9cSBarry Smith . pc - the multigrid context 8293fc8bf9cSBarry Smith 8303fc8bf9cSBarry Smith Output Parameter: 8313fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 8323fc8bf9cSBarry Smith 8333fc8bf9cSBarry Smith Options Database Key: 8343fc8bf9cSBarry Smith $ -pc_mg_galerkin 8353fc8bf9cSBarry Smith 8363fc8bf9cSBarry Smith Level: intermediate 8373fc8bf9cSBarry Smith 8383fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 8393fc8bf9cSBarry Smith 8403fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 8413fc8bf9cSBarry Smith 8423fc8bf9cSBarry Smith @*/ 8433fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 8443fc8bf9cSBarry Smith { 845f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 846f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 8473fc8bf9cSBarry Smith 8483fc8bf9cSBarry Smith PetscFunctionBegin; 8493fc8bf9cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 850f3fbd535SBarry Smith if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 851f3fbd535SBarry Smith *galerkin = mglevels[0]->galerkin; 8523fc8bf9cSBarry Smith PetscFunctionReturn(0); 8533fc8bf9cSBarry Smith } 8543fc8bf9cSBarry Smith 8553fc8bf9cSBarry Smith #undef __FUNCT__ 8569dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 8574b9ad928SBarry Smith /*@ 85897177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 85997177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 8604b9ad928SBarry Smith pre-smoothing steps on different levels. 8614b9ad928SBarry Smith 8624b9ad928SBarry Smith Collective on PC 8634b9ad928SBarry Smith 8644b9ad928SBarry Smith Input Parameters: 8654b9ad928SBarry Smith + mg - the multigrid context 8664b9ad928SBarry Smith - n - the number of smoothing steps 8674b9ad928SBarry Smith 8684b9ad928SBarry Smith Options Database Key: 8694b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 8704b9ad928SBarry Smith 8714b9ad928SBarry Smith Level: advanced 8724b9ad928SBarry Smith 8734b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 8744b9ad928SBarry Smith 87597177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 8764b9ad928SBarry Smith @*/ 87797177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n) 8784b9ad928SBarry Smith { 879f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 880f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 8816849ba73SBarry Smith PetscErrorCode ierr; 88279416396SBarry Smith PetscInt i,levels; 8834b9ad928SBarry Smith 8844b9ad928SBarry Smith PetscFunctionBegin; 8854482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 886f3fbd535SBarry Smith if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 887f3fbd535SBarry Smith levels = mglevels[0]->levels; 8884b9ad928SBarry Smith 889b05257ddSBarry Smith for (i=1; i<levels; i++) { 8904b9ad928SBarry Smith /* make sure smoother up and down are different */ 89197177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 892f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 89331567311SBarry Smith mg->default_smoothd = n; 8944b9ad928SBarry Smith } 8954b9ad928SBarry Smith PetscFunctionReturn(0); 8964b9ad928SBarry Smith } 8974b9ad928SBarry Smith 8984b9ad928SBarry Smith #undef __FUNCT__ 8999dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 9004b9ad928SBarry Smith /*@ 90197177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 90297177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 9034b9ad928SBarry Smith post-smoothing steps on different levels. 9044b9ad928SBarry Smith 9054b9ad928SBarry Smith Collective on PC 9064b9ad928SBarry Smith 9074b9ad928SBarry Smith Input Parameters: 9084b9ad928SBarry Smith + mg - the multigrid context 9094b9ad928SBarry Smith - n - the number of smoothing steps 9104b9ad928SBarry Smith 9114b9ad928SBarry Smith Options Database Key: 9124b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 9134b9ad928SBarry Smith 9144b9ad928SBarry Smith Level: advanced 9154b9ad928SBarry Smith 9164b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 917a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 9184b9ad928SBarry Smith 9194b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 9204b9ad928SBarry Smith 92197177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 9224b9ad928SBarry Smith @*/ 92397177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 9244b9ad928SBarry Smith { 925f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 926f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 9276849ba73SBarry Smith PetscErrorCode ierr; 92879416396SBarry Smith PetscInt i,levels; 9294b9ad928SBarry Smith 9304b9ad928SBarry Smith PetscFunctionBegin; 9314482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 932f3fbd535SBarry Smith if (!mglevels) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 933f3fbd535SBarry Smith levels = mglevels[0]->levels; 9344b9ad928SBarry Smith 9354b9ad928SBarry Smith for (i=1; i<levels; i++) { 9364b9ad928SBarry Smith /* make sure smoother up and down are different */ 93797177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 938f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 93931567311SBarry Smith mg->default_smoothu = n; 9404b9ad928SBarry Smith } 9414b9ad928SBarry Smith PetscFunctionReturn(0); 9424b9ad928SBarry Smith } 9434b9ad928SBarry Smith 9444b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 9454b9ad928SBarry Smith 9463b09bd56SBarry Smith /*MC 947ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 9483b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 9493b09bd56SBarry Smith 9503b09bd56SBarry Smith Options Database Keys: 9513b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 9520d353602SBarry Smith . -pc_mg_cycles v or w 95379416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 9543b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 9553b09bd56SBarry Smith . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 9563b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 9573b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 95868eff7e6SBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators 9593b09bd56SBarry Smith - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 9603b09bd56SBarry Smith to the Socket viewer for reading from Matlab. 9613b09bd56SBarry Smith 9623b09bd56SBarry Smith Notes: 9633b09bd56SBarry Smith 9643b09bd56SBarry Smith Level: intermediate 9653b09bd56SBarry Smith 9668f87f92bSBarry Smith Concepts: multigrid/multilevel 9673b09bd56SBarry Smith 9683b09bd56SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 9690d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 97097177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 97197177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 9720d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 9733b09bd56SBarry Smith M*/ 9743b09bd56SBarry Smith 9754b9ad928SBarry Smith EXTERN_C_BEGIN 9764b9ad928SBarry Smith #undef __FUNCT__ 9774b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 978dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 9794b9ad928SBarry Smith { 980f3fbd535SBarry Smith PC_MG *mg; 981f3fbd535SBarry Smith PetscErrorCode ierr; 982f3fbd535SBarry Smith 9834b9ad928SBarry Smith PetscFunctionBegin; 984f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 985f3fbd535SBarry Smith pc->data = (void*)mg; 986f3fbd535SBarry Smith mg->nlevels = -1; 987f3fbd535SBarry Smith 9884b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 9894b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 9904b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 9914b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 9924b9ad928SBarry Smith pc->ops->view = PCView_MG; 9934b9ad928SBarry Smith PetscFunctionReturn(0); 9944b9ad928SBarry Smith } 9954b9ad928SBarry Smith EXTERN_C_END 996