1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines the multigrid preconditioner interface. 44b9ad928SBarry Smith */ 5*c6db04a5SJed Brown #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscpcmg.h" I*/ 64b9ad928SBarry Smith 74b9ad928SBarry Smith 84b9ad928SBarry Smith #undef __FUNCT__ 99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private" 1031567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 114b9ad928SBarry Smith { 1231567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1331567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 146849ba73SBarry Smith PetscErrorCode ierr; 1531567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 164b9ad928SBarry Smith 174b9ad928SBarry Smith PetscFunctionBegin; 184b9ad928SBarry Smith 1963e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2031567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 2163e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2231567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 2363e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 2431567311SBarry Smith ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 2563e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 264b9ad928SBarry Smith 274b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 2831567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 294b9ad928SBarry Smith PetscReal rnorm; 3031567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 314b9ad928SBarry Smith if (rnorm <= mg->ttol) { 3270441072SBarry Smith if (rnorm < mg->abstol) { 334d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 341e2582c4SBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %G is less than absolute tolerance %G\n",rnorm,mg->abstol);CHKERRQ(ierr); 354b9ad928SBarry Smith } else { 364d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 371e2582c4SBarry 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); 384b9ad928SBarry Smith } 394b9ad928SBarry Smith PetscFunctionReturn(0); 404b9ad928SBarry Smith } 414b9ad928SBarry Smith } 424b9ad928SBarry Smith 4331567311SBarry Smith mgc = *(mglevelsin - 1); 4463e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 4531567311SBarry Smith ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 4663e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 47efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 484b9ad928SBarry Smith while (cycles--) { 4931567311SBarry Smith ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 504b9ad928SBarry Smith } 5163e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5231567311SBarry Smith ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 5363e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5463e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 5531567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 5663e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 574b9ad928SBarry Smith } 584b9ad928SBarry Smith PetscFunctionReturn(0); 594b9ad928SBarry Smith } 604b9ad928SBarry Smith 614b9ad928SBarry Smith #undef __FUNCT__ 624b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG" 63ace3abfcSBarry Smith 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) 644b9ad928SBarry Smith { 65f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 66f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 67dfbe8321SBarry Smith PetscErrorCode ierr; 68f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 694b9ad928SBarry Smith 704b9ad928SBarry Smith PetscFunctionBegin; 71f3fbd535SBarry Smith mglevels[levels-1]->b = b; 72f3fbd535SBarry Smith mglevels[levels-1]->x = x; 734b9ad928SBarry Smith 7431567311SBarry Smith mg->rtol = rtol; 7531567311SBarry Smith mg->abstol = abstol; 7631567311SBarry Smith mg->dtol = dtol; 774b9ad928SBarry Smith if (rtol) { 784b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 794b9ad928SBarry Smith PetscReal rnorm; 807319c654SBarry Smith if (zeroguess) { 817319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 827319c654SBarry Smith } else { 83f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 844b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 857319c654SBarry Smith } 8631567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 8770441072SBarry Smith } else if (abstol) { 8831567311SBarry Smith mg->ttol = abstol; 894b9ad928SBarry Smith } else { 9031567311SBarry Smith mg->ttol = 0.0; 914b9ad928SBarry Smith } 924b9ad928SBarry Smith 934d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 944d0a8057SBarry Smith stop prematurely do to small residual */ 954d0a8057SBarry Smith for (i=1; i<levels; i++) { 96f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 97f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 98f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 994b9ad928SBarry Smith } 1004d0a8057SBarry Smith } 1014d0a8057SBarry Smith 1024d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1034d0a8057SBarry Smith for (i=0; i<its; i++) { 104f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1054d0a8057SBarry Smith if (*reason) break; 1064d0a8057SBarry Smith } 1074d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1084d0a8057SBarry Smith *outits = i; 1094b9ad928SBarry Smith PetscFunctionReturn(0); 1104b9ad928SBarry Smith } 1114b9ad928SBarry Smith 1124b9ad928SBarry Smith #undef __FUNCT__ 1139dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 1144b9ad928SBarry Smith /*@C 11597177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 1164b9ad928SBarry Smith Must be called before any other MG routine. 1174b9ad928SBarry Smith 118ad4df100SBarry Smith Logically Collective on PC 1194b9ad928SBarry Smith 1204b9ad928SBarry Smith Input Parameters: 1214b9ad928SBarry Smith + pc - the preconditioner context 1224b9ad928SBarry Smith . levels - the number of levels 1234b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 1244b9ad928SBarry Smith on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 1254b9ad928SBarry Smith 1264b9ad928SBarry Smith Level: intermediate 1274b9ad928SBarry Smith 1284b9ad928SBarry Smith Notes: 1294b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 1304b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 1314b9ad928SBarry Smith 1324b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 1334b9ad928SBarry Smith 13497177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 1354b9ad928SBarry Smith @*/ 1367087cfbeSBarry Smith PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 1374b9ad928SBarry Smith { 138dfbe8321SBarry Smith PetscErrorCode ierr; 139f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 140f3fbd535SBarry Smith MPI_Comm comm = ((PetscObject)pc)->comm; 141f3fbd535SBarry Smith PC_MG_Levels **mglevels; 142f3fbd535SBarry Smith PetscInt i; 143f3fbd535SBarry Smith PetscMPIInt size; 144f3fbd535SBarry Smith const char *prefix; 145f3fbd535SBarry Smith PC ipc; 1464b9ad928SBarry Smith 1474b9ad928SBarry Smith PetscFunctionBegin; 1480700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 149e7e72b3dSBarry Smith if (mg->nlevels > -1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Number levels already set for MG\n make sure that you call PCMGSetLevels() before KSPSetFromOptions()"); 150e32f2f54SBarry Smith if (mg->levels) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Internal error in PETSc, this array should not yet exist"); 151c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,levels,2); 152f3fbd535SBarry Smith 153f3fbd535SBarry Smith mg->nlevels = levels; 154218a07d4SBarry Smith mg->galerkin = PETSC_FALSE; 155218a07d4SBarry Smith mg->galerkinused = PETSC_FALSE; 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; 16731567311SBarry Smith mg->default_smoothu = 1; 16831567311SBarry Smith mg->default_smoothd = 1; 16963e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 17063e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 17163e6d426SJed Brown mglevels[i]->eventresidual = 0; 17263e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 173f3fbd535SBarry Smith 174f3fbd535SBarry Smith if (comms) comm = comms[i]; 175f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 176f3fbd535SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 17731567311SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 178f3fbd535SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 179f3fbd535SBarry Smith 180f3fbd535SBarry Smith /* do special stuff for coarse grid */ 181f3fbd535SBarry Smith if (!i && levels > 1) { 182f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 183f3fbd535SBarry Smith 184f3fbd535SBarry Smith /* coarse solve is (redundant) LU by default */ 185f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 186f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 187f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 188f3fbd535SBarry Smith if (size > 1) { 189f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 190f3fbd535SBarry Smith } else { 191f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 192f3fbd535SBarry Smith } 193f3fbd535SBarry Smith 194f3fbd535SBarry Smith } else { 195f3fbd535SBarry Smith char tprefix[128]; 196f3fbd535SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 197f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 198f3fbd535SBarry Smith } 199f3fbd535SBarry Smith ierr = PetscLogObjectParent(pc,mglevels[i]->smoothd);CHKERRQ(ierr); 200f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 20131567311SBarry Smith mg->rtol = 0.0; 20231567311SBarry Smith mg->abstol = 0.0; 20331567311SBarry Smith mg->dtol = 0.0; 20431567311SBarry Smith mg->ttol = 0.0; 20531567311SBarry Smith mg->cyclesperpcapply = 1; 206f3fbd535SBarry Smith } 20731567311SBarry Smith mg->am = PC_MG_MULTIPLICATIVE; 208f3fbd535SBarry Smith mg->levels = mglevels; 2094b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_MG; 2104b9ad928SBarry Smith PetscFunctionReturn(0); 2114b9ad928SBarry Smith } 2124b9ad928SBarry Smith 2134b9ad928SBarry Smith #undef __FUNCT__ 214c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG_Private" 215c07bf074SBarry Smith PetscErrorCode PCDestroy_MG_Private(PC pc) 216f3fbd535SBarry Smith { 217f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 218f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 219f3fbd535SBarry Smith PetscErrorCode ierr; 220f3fbd535SBarry Smith PetscInt i,n; 221f3fbd535SBarry Smith 222f3fbd535SBarry Smith PetscFunctionBegin; 223f3fbd535SBarry Smith if (mglevels) { 224f3fbd535SBarry Smith n = mglevels[0]->levels; 225f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 226f3fbd535SBarry Smith if (mglevels[i+1]->r) {ierr = VecDestroy(mglevels[i+1]->r);CHKERRQ(ierr);} 227f3fbd535SBarry Smith if (mglevels[i]->b) {ierr = VecDestroy(mglevels[i]->b);CHKERRQ(ierr);} 228f3fbd535SBarry Smith if (mglevels[i]->x) {ierr = VecDestroy(mglevels[i]->x);CHKERRQ(ierr);} 229f3fbd535SBarry Smith if (mglevels[i+1]->restrct) {ierr = MatDestroy(mglevels[i+1]->restrct);CHKERRQ(ierr);} 230f3fbd535SBarry Smith if (mglevels[i+1]->interpolate) {ierr = MatDestroy(mglevels[i+1]->interpolate);CHKERRQ(ierr);} 231f3fbd535SBarry Smith } 232f3fbd535SBarry Smith 233f3fbd535SBarry Smith for (i=0; i<n; i++) { 234f3fbd535SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 235f3fbd535SBarry Smith ierr = KSPDestroy(mglevels[i]->smoothd);CHKERRQ(ierr); 236f3fbd535SBarry Smith } 237f3fbd535SBarry Smith ierr = KSPDestroy(mglevels[i]->smoothu);CHKERRQ(ierr); 238f3fbd535SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 239f3fbd535SBarry Smith } 240c31cb41cSBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 241f3fbd535SBarry Smith } 242c07bf074SBarry Smith mg->nlevels = -1; 243c07bf074SBarry Smith PetscFunctionReturn(0); 244c07bf074SBarry Smith } 245c07bf074SBarry Smith 246c07bf074SBarry Smith #undef __FUNCT__ 247c07bf074SBarry Smith #define __FUNCT__ "PCDestroy_MG" 248c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 249c07bf074SBarry Smith { 250c07bf074SBarry Smith PetscErrorCode ierr; 251c07bf074SBarry Smith 252c07bf074SBarry Smith PetscFunctionBegin; 253c07bf074SBarry Smith ierr = PCDestroy_MG_Private(pc);CHKERRQ(ierr); 254c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 255f3fbd535SBarry Smith PetscFunctionReturn(0); 256f3fbd535SBarry Smith } 257f3fbd535SBarry Smith 258f3fbd535SBarry Smith 259f3fbd535SBarry Smith 26009573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 26109573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 26209573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 263f3fbd535SBarry Smith 264f3fbd535SBarry Smith /* 265f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 266f3fbd535SBarry Smith or full cycle of multigrid. 267f3fbd535SBarry Smith 268f3fbd535SBarry Smith Note: 269f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 270f3fbd535SBarry Smith */ 271f3fbd535SBarry Smith #undef __FUNCT__ 272f3fbd535SBarry Smith #define __FUNCT__ "PCApply_MG" 273f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 274f3fbd535SBarry Smith { 275f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 276f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 277f3fbd535SBarry Smith PetscErrorCode ierr; 278f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 279f3fbd535SBarry Smith 280f3fbd535SBarry Smith PetscFunctionBegin; 281e1d8e5deSBarry Smith 282e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 283e1d8e5deSBarry Smith for (i=0; i<levels-1; i++) { 284e1d8e5deSBarry Smith if (!mglevels[i]->A) { 285e1d8e5deSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 286e1d8e5deSBarry Smith } 287e1d8e5deSBarry Smith } 288e1d8e5deSBarry Smith 289f3fbd535SBarry Smith mglevels[levels-1]->b = b; 290f3fbd535SBarry Smith mglevels[levels-1]->x = x; 29131567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 292f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 29331567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 294f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_NULL);CHKERRQ(ierr); 295f3fbd535SBarry Smith } 296f3fbd535SBarry Smith } 29731567311SBarry Smith else if (mg->am == PC_MG_ADDITIVE) { 29831567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 299f3fbd535SBarry Smith } 30031567311SBarry Smith else if (mg->am == PC_MG_KASKADE) { 30131567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 302f3fbd535SBarry Smith } 303f3fbd535SBarry Smith else { 304f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 305f3fbd535SBarry Smith } 306f3fbd535SBarry Smith PetscFunctionReturn(0); 307f3fbd535SBarry Smith } 308f3fbd535SBarry Smith 309f3fbd535SBarry Smith 310f3fbd535SBarry Smith #undef __FUNCT__ 311f3fbd535SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 312f3fbd535SBarry Smith PetscErrorCode PCSetFromOptions_MG(PC pc) 313f3fbd535SBarry Smith { 314f3fbd535SBarry Smith PetscErrorCode ierr; 315f3fbd535SBarry Smith PetscInt m,levels = 1,cycles; 316ace3abfcSBarry Smith PetscBool flg; 317f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 318f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 319f3fbd535SBarry Smith PCMGType mgtype; 320f3fbd535SBarry Smith PCMGCycleType mgctype; 321f3fbd535SBarry Smith 322f3fbd535SBarry Smith PetscFunctionBegin; 323f3fbd535SBarry Smith ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 32418aabeadSBarry Smith if (!mglevels) { 325f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 326f3fbd535SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 327f3fbd535SBarry Smith mglevels = mg->levels; 328f3fbd535SBarry Smith } 329f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 330f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 331f3fbd535SBarry Smith if (flg) { 332f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 333f3fbd535SBarry Smith }; 334f3fbd535SBarry Smith flg = PETSC_FALSE; 335acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 336f3fbd535SBarry Smith if (flg) { 337f3fbd535SBarry Smith ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 338f3fbd535SBarry Smith } 339f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 340f3fbd535SBarry Smith if (flg) { 341f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 342f3fbd535SBarry Smith } 343f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 344f3fbd535SBarry Smith if (flg) { 345f3fbd535SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 346f3fbd535SBarry Smith } 34731567311SBarry Smith mgtype = mg->am; 348f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 349f3fbd535SBarry Smith if (flg) { 350f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 351f3fbd535SBarry Smith } 35231567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 35331567311SBarry Smith ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 354f3fbd535SBarry Smith if (flg) { 355f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 356f3fbd535SBarry Smith } 357f3fbd535SBarry Smith } 358f3fbd535SBarry Smith flg = PETSC_FALSE; 359acfcf0e5SJed Brown ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,PETSC_NULL);CHKERRQ(ierr); 360f3fbd535SBarry Smith if (flg) { 361f3fbd535SBarry Smith PetscInt i; 362f3fbd535SBarry Smith char eventname[128]; 36363e6d426SJed Brown if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 364f3fbd535SBarry Smith levels = mglevels[0]->levels; 365f3fbd535SBarry Smith for (i=0; i<levels; i++) { 366f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 36763e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 368f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 36963e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 370f3fbd535SBarry Smith if (i) { 371f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 37263e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 373f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 37463e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 375f3fbd535SBarry Smith } 376f3fbd535SBarry Smith } 377f3fbd535SBarry Smith } 378f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 379f3fbd535SBarry Smith PetscFunctionReturn(0); 380f3fbd535SBarry Smith } 381f3fbd535SBarry Smith 382f3fbd535SBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 383f3fbd535SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 384f3fbd535SBarry Smith 385f3fbd535SBarry Smith #undef __FUNCT__ 386f3fbd535SBarry Smith #define __FUNCT__ "PCView_MG" 387f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 388f3fbd535SBarry Smith { 389f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 390f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 391f3fbd535SBarry Smith PetscErrorCode ierr; 392f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 393ace3abfcSBarry Smith PetscBool iascii; 394f3fbd535SBarry Smith 395f3fbd535SBarry Smith PetscFunctionBegin; 3962692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 397f3fbd535SBarry Smith if (iascii) { 39831567311SBarry 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); 39931567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 40031567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 401f3fbd535SBarry Smith } 402218a07d4SBarry Smith if (mg->galerkin) { 403f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 4044f66f45eSBarry Smith } else { 4054f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 406f3fbd535SBarry Smith } 407f3fbd535SBarry Smith for (i=0; i<levels; i++) { 408f3fbd535SBarry Smith if (!i) { 4092d19a89fSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr); 410f3fbd535SBarry Smith } else { 41131567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothd);CHKERRQ(ierr); 412f3fbd535SBarry Smith } 413f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 414f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 415f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 416f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 417f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 418f3fbd535SBarry Smith } else if (i){ 41931567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg->default_smoothu);CHKERRQ(ierr); 420f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 421f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 422f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 423f3fbd535SBarry Smith } 424f3fbd535SBarry Smith } 425f3fbd535SBarry Smith } else { 42665e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 427f3fbd535SBarry Smith } 428f3fbd535SBarry Smith PetscFunctionReturn(0); 429f3fbd535SBarry Smith } 430f3fbd535SBarry Smith 431f3fbd535SBarry Smith /* 432f3fbd535SBarry Smith Calls setup for the KSP on each level 433f3fbd535SBarry Smith */ 434f3fbd535SBarry Smith #undef __FUNCT__ 435f3fbd535SBarry Smith #define __FUNCT__ "PCSetUp_MG" 436f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 437f3fbd535SBarry Smith { 438f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 439f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 440f3fbd535SBarry Smith PetscErrorCode ierr; 441f3fbd535SBarry Smith PetscInt i,n = mglevels[0]->levels; 442f3fbd535SBarry Smith PC cpc,mpc; 443ace3abfcSBarry Smith PetscBool preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump = PETSC_FALSE,opsset; 444f3fbd535SBarry Smith PetscViewerASCIIMonitor ascii; 445f3fbd535SBarry Smith PetscViewer viewer = PETSC_NULL; 446f3fbd535SBarry Smith MPI_Comm comm; 447f3fbd535SBarry Smith Mat dA,dB; 448f3fbd535SBarry Smith MatStructure uflag; 449f3fbd535SBarry Smith Vec tvec; 450218a07d4SBarry Smith DM *dms; 451f3fbd535SBarry Smith 452f3fbd535SBarry Smith PetscFunctionBegin; 453f3fbd535SBarry Smith 454f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 455f3fbd535SBarry Smith /* so use those from global PC */ 456f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 457f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 458f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 459f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[n-1]->smoothd,&mpc);CHKERRQ(ierr); 4601338a6b9SJed Brown if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2)) || ((mpc == cpc) && (mpc->setupcalled == 2))) { 461f3fbd535SBarry 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); 462f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 463f3fbd535SBarry Smith } 464f3fbd535SBarry Smith 4652d2b81a6SBarry Smith if (pc->dm && !pc->setupcalled) { 4662d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 4672e499ae9SBarry Smith Mat p; 468218a07d4SBarry Smith ierr = PetscMalloc(n*sizeof(DM),&dms);CHKERRQ(ierr); 469218a07d4SBarry Smith dms[n-1] = pc->dm; 470218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 471218a07d4SBarry Smith ierr = DMCoarsen(dms[i+1],PETSC_NULL,&dms[i]);CHKERRQ(ierr); 47211629dbeSBarry Smith ierr = DMSetFunction(dms[i],0); 47311629dbeSBarry Smith ierr = DMSetInitialGuess(dms[i],0); 47424c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 4752d2b81a6SBarry Smith ierr = DMGetInterpolation(dms[i],dms[i+1],&p,PETSC_NULL);CHKERRQ(ierr); 4762d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 4772d2b81a6SBarry Smith ierr = MatDestroy(p);CHKERRQ(ierr); 478218a07d4SBarry Smith } 47924c3aa18SBarry Smith } 4802d2b81a6SBarry Smith 4812d2b81a6SBarry Smith if (!mg->galerkin) { 482d42688cbSBarry Smith /* each coarse level gets its DM; finest level does not get DM because it shared the outer PC operators */ 4832d2b81a6SBarry Smith for (i=n-2; i>-1; i--) { 4842e499ae9SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 4852d2b81a6SBarry Smith } 4862d2b81a6SBarry Smith } 4872d2b81a6SBarry Smith 488218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 489218a07d4SBarry Smith ierr = DMDestroy(dms[i]);CHKERRQ(ierr); 490218a07d4SBarry Smith } 4912d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 492218a07d4SBarry Smith } 493218a07d4SBarry Smith 494218a07d4SBarry Smith if (mg->galerkin) { 495f3fbd535SBarry Smith Mat B; 496218a07d4SBarry Smith mg->galerkinused = PETSC_TRUE; 497f3fbd535SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 498f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 499f3fbd535SBarry Smith if (!pc->setupcalled) { 500f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 501f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 502f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 503f3fbd535SBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 504f3fbd535SBarry Smith dB = B; 505f3fbd535SBarry Smith } 506cd9507b2SBarry Smith if (n > 1) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 507f3fbd535SBarry Smith } else { 508f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 509f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 510f3fbd535SBarry Smith ierr = MatPtAP(dB,mglevels[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 511f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 512f3fbd535SBarry Smith dB = B; 513f3fbd535SBarry Smith } 514f3fbd535SBarry Smith } 515f3fbd535SBarry Smith } 516f3fbd535SBarry Smith 517f3fbd535SBarry Smith if (!pc->setupcalled) { 518acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_monitor",&monitor,PETSC_NULL);CHKERRQ(ierr); 519f3fbd535SBarry Smith 520f3fbd535SBarry Smith for (i=0; i<n; i++) { 521f3fbd535SBarry Smith if (monitor) { 522f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothd,&comm);CHKERRQ(ierr); 523f3fbd535SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 524f3fbd535SBarry Smith ierr = KSPMonitorSet(mglevels[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 525f3fbd535SBarry Smith } 526f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 527f3fbd535SBarry Smith } 528f3fbd535SBarry Smith for (i=1; i<n; i++) { 529f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 530f3fbd535SBarry Smith if (monitor) { 531f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[i]->smoothu,&comm);CHKERRQ(ierr); 532f3fbd535SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 533f3fbd535SBarry Smith ierr = KSPMonitorSet(mglevels[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 534f3fbd535SBarry Smith } 535f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 536f3fbd535SBarry Smith } 537f3fbd535SBarry Smith } 538f3fbd535SBarry Smith for (i=1; i<n; i++) { 539f3fbd535SBarry Smith if (mglevels[i]->restrct && !mglevels[i]->interpolate) { 540f3fbd535SBarry Smith ierr = PCMGSetInterpolation(pc,i,mglevels[i]->restrct);CHKERRQ(ierr); 541f3fbd535SBarry Smith } 542f3fbd535SBarry Smith if (!mglevels[i]->restrct && mglevels[i]->interpolate) { 543f3fbd535SBarry Smith ierr = PCMGSetRestriction(pc,i,mglevels[i]->interpolate);CHKERRQ(ierr); 544f3fbd535SBarry Smith } 545f3fbd535SBarry Smith #if defined(PETSC_USE_DEBUG) 546f3fbd535SBarry Smith if (!mglevels[i]->restrct || !mglevels[i]->interpolate) { 54765e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 548f3fbd535SBarry Smith } 549f3fbd535SBarry Smith #endif 550f3fbd535SBarry Smith } 551f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 552f3fbd535SBarry Smith if (!mglevels[i]->b) { 553f3fbd535SBarry Smith Vec *vec; 554f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 555f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 556f3fbd535SBarry Smith ierr = VecDestroy(*vec);CHKERRQ(ierr); 557f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 558f3fbd535SBarry Smith } 559f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 560f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 561f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 562f3fbd535SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 563f3fbd535SBarry Smith } 564f3fbd535SBarry Smith if (!mglevels[i]->x) { 565f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 566f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 567f3fbd535SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 568f3fbd535SBarry Smith } 569f3fbd535SBarry Smith } 570f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 571f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 572f3fbd535SBarry Smith Vec *vec; 573f3fbd535SBarry Smith ierr = KSPGetVecs(mglevels[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 574f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 575f3fbd535SBarry Smith ierr = VecDestroy(*vec);CHKERRQ(ierr); 576f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 577f3fbd535SBarry Smith } 578f3fbd535SBarry Smith } 579f3fbd535SBarry Smith 580f3fbd535SBarry Smith 581f3fbd535SBarry Smith for (i=1; i<n; i++) { 582f3fbd535SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd) { 583f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 584f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 585f3fbd535SBarry Smith } 58663e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 587f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 58863e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 589d42688cbSBarry Smith if (!mglevels[i]->residual) { 590d42688cbSBarry Smith Mat mat; 591d42688cbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 592d42688cbSBarry Smith ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 593d42688cbSBarry Smith } 594f3fbd535SBarry Smith } 595f3fbd535SBarry Smith for (i=1; i<n; i++) { 596f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 597f3fbd535SBarry Smith Mat downmat,downpmat; 598f3fbd535SBarry Smith MatStructure matflag; 599ace3abfcSBarry Smith PetscBool opsset; 600f3fbd535SBarry Smith 601f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 602f3fbd535SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 603f3fbd535SBarry Smith if (!opsset) { 604f3fbd535SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 605f3fbd535SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 606f3fbd535SBarry Smith } 607f3fbd535SBarry Smith 608f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 60963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 610f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 61163e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 612f3fbd535SBarry Smith } 613f3fbd535SBarry Smith } 614f3fbd535SBarry Smith 615f3fbd535SBarry Smith /* 616f3fbd535SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 617f3fbd535SBarry Smith */ 618f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)mglevels[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 619f3fbd535SBarry Smith if (preonly) { 620f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 621f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 622f3fbd535SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 623f3fbd535SBarry Smith if (!lu && !redundant && !cholesky) { 624f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 625f3fbd535SBarry Smith } 626f3fbd535SBarry Smith } 627f3fbd535SBarry Smith 628f3fbd535SBarry Smith if (!pc->setupcalled) { 629f3fbd535SBarry Smith if (monitor) { 630f3fbd535SBarry Smith ierr = PetscObjectGetComm((PetscObject)mglevels[0]->smoothd,&comm);CHKERRQ(ierr); 631f3fbd535SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 632f3fbd535SBarry Smith ierr = KSPMonitorSet(mglevels[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 633f3fbd535SBarry Smith } 634f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[0]->smoothd);CHKERRQ(ierr); 635f3fbd535SBarry Smith } 636f3fbd535SBarry Smith 63763e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 638f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 63963e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 640f3fbd535SBarry Smith 641f3fbd535SBarry Smith /* 642f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 643e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 644f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 645f3fbd535SBarry Smith 646f3fbd535SBarry Smith Only support one or the other at the same time. 647f3fbd535SBarry Smith */ 648f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 649acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,PETSC_NULL);CHKERRQ(ierr); 650f3fbd535SBarry Smith if (dump) { 651f3fbd535SBarry Smith viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 652f3fbd535SBarry Smith } 653f3fbd535SBarry Smith dump = PETSC_FALSE; 654f3fbd535SBarry Smith #endif 655acfcf0e5SJed Brown ierr = PetscOptionsGetBool(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,PETSC_NULL);CHKERRQ(ierr); 656f3fbd535SBarry Smith if (dump) { 657f3fbd535SBarry Smith viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 658f3fbd535SBarry Smith } 659f3fbd535SBarry Smith 660f3fbd535SBarry Smith if (viewer) { 661f3fbd535SBarry Smith for (i=1; i<n; i++) { 662f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 663f3fbd535SBarry Smith } 664f3fbd535SBarry Smith for (i=0; i<n; i++) { 665f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 666f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 667f3fbd535SBarry Smith } 668f3fbd535SBarry Smith } 669f3fbd535SBarry Smith PetscFunctionReturn(0); 670f3fbd535SBarry Smith } 671f3fbd535SBarry Smith 672f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 673f3fbd535SBarry Smith 674f3fbd535SBarry Smith #undef __FUNCT__ 6759dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 6764b9ad928SBarry Smith /*@ 67797177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 6784b9ad928SBarry Smith 6794b9ad928SBarry Smith Not Collective 6804b9ad928SBarry Smith 6814b9ad928SBarry Smith Input Parameter: 6824b9ad928SBarry Smith . pc - the preconditioner context 6834b9ad928SBarry Smith 6844b9ad928SBarry Smith Output parameter: 6854b9ad928SBarry Smith . levels - the number of levels 6864b9ad928SBarry Smith 6874b9ad928SBarry Smith Level: advanced 6884b9ad928SBarry Smith 6894b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 6904b9ad928SBarry Smith 69197177400SBarry Smith .seealso: PCMGSetLevels() 6924b9ad928SBarry Smith @*/ 6937087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 6944b9ad928SBarry Smith { 695f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 6964b9ad928SBarry Smith 6974b9ad928SBarry Smith PetscFunctionBegin; 6980700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 6994482741eSBarry Smith PetscValidIntPointer(levels,2); 700f3fbd535SBarry Smith *levels = mg->nlevels; 7014b9ad928SBarry Smith PetscFunctionReturn(0); 7024b9ad928SBarry Smith } 7034b9ad928SBarry Smith 7044b9ad928SBarry Smith #undef __FUNCT__ 7059dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 7064b9ad928SBarry Smith /*@ 70797177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 7084b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 7094b9ad928SBarry Smith 710ad4df100SBarry Smith Logically Collective on PC 7114b9ad928SBarry Smith 7124b9ad928SBarry Smith Input Parameters: 7134b9ad928SBarry Smith + pc - the preconditioner context 7149dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 7159dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 7164b9ad928SBarry Smith 7174b9ad928SBarry Smith Options Database Key: 7184b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 7194b9ad928SBarry Smith additive, full, kaskade 7204b9ad928SBarry Smith 7214b9ad928SBarry Smith Level: advanced 7224b9ad928SBarry Smith 7234b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 7244b9ad928SBarry Smith 72597177400SBarry Smith .seealso: PCMGSetLevels() 7264b9ad928SBarry Smith @*/ 7277087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 7284b9ad928SBarry Smith { 729f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 7304b9ad928SBarry Smith 7314b9ad928SBarry Smith PetscFunctionBegin; 7320700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 733c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 73431567311SBarry Smith mg->am = form; 7359dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 7364b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 7374b9ad928SBarry Smith PetscFunctionReturn(0); 7384b9ad928SBarry Smith } 7394b9ad928SBarry Smith 7404b9ad928SBarry Smith #undef __FUNCT__ 7410d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 7424b9ad928SBarry Smith /*@ 7430d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 7444b9ad928SBarry Smith complicated cycling. 7454b9ad928SBarry Smith 746ad4df100SBarry Smith Logically Collective on PC 7474b9ad928SBarry Smith 7484b9ad928SBarry Smith Input Parameters: 749c2be2410SBarry Smith + pc - the multigrid context 7500d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 7514b9ad928SBarry Smith 7524b9ad928SBarry Smith Options Database Key: 7530d353602SBarry Smith $ -pc_mg_cycle_type v or w 7544b9ad928SBarry Smith 7554b9ad928SBarry Smith Level: advanced 7564b9ad928SBarry Smith 7574b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 7584b9ad928SBarry Smith 7590d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 7604b9ad928SBarry Smith @*/ 7617087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 7624b9ad928SBarry Smith { 763f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 764f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 76579416396SBarry Smith PetscInt i,levels; 7664b9ad928SBarry Smith 7674b9ad928SBarry Smith PetscFunctionBegin; 7680700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 769e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 770c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 771f3fbd535SBarry Smith levels = mglevels[0]->levels; 7724b9ad928SBarry Smith 7734b9ad928SBarry Smith for (i=0; i<levels; i++) { 774f3fbd535SBarry Smith mglevels[i]->cycles = n; 7754b9ad928SBarry Smith } 7764b9ad928SBarry Smith PetscFunctionReturn(0); 7774b9ad928SBarry Smith } 7784b9ad928SBarry Smith 7794b9ad928SBarry Smith #undef __FUNCT__ 7808cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 7818cc2d5dfSBarry Smith /*@ 7828cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 7838cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 7848cc2d5dfSBarry Smith 785ad4df100SBarry Smith Logically Collective on PC 7868cc2d5dfSBarry Smith 7878cc2d5dfSBarry Smith Input Parameters: 7888cc2d5dfSBarry Smith + pc - the multigrid context 7898cc2d5dfSBarry Smith - n - number of cycles (default is 1) 7908cc2d5dfSBarry Smith 7918cc2d5dfSBarry Smith Options Database Key: 7928cc2d5dfSBarry Smith $ -pc_mg_multiplicative_cycles n 7938cc2d5dfSBarry Smith 7948cc2d5dfSBarry Smith Level: advanced 7958cc2d5dfSBarry Smith 7968cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 7978cc2d5dfSBarry Smith 7988cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 7998cc2d5dfSBarry Smith 8008cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 8018cc2d5dfSBarry Smith @*/ 8027087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 8038cc2d5dfSBarry Smith { 804f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 805f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 8068cc2d5dfSBarry Smith PetscInt i,levels; 8078cc2d5dfSBarry Smith 8088cc2d5dfSBarry Smith PetscFunctionBegin; 8090700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 810e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 811c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 812f3fbd535SBarry Smith levels = mglevels[0]->levels; 8138cc2d5dfSBarry Smith 8148cc2d5dfSBarry Smith for (i=0; i<levels; i++) { 81531567311SBarry Smith mg->cyclesperpcapply = n; 8168cc2d5dfSBarry Smith } 8178cc2d5dfSBarry Smith PetscFunctionReturn(0); 8188cc2d5dfSBarry Smith } 8198cc2d5dfSBarry Smith 8208cc2d5dfSBarry Smith #undef __FUNCT__ 8219dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 822c2be2410SBarry Smith /*@ 82397177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 824c2be2410SBarry Smith finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 825c2be2410SBarry Smith 826ad4df100SBarry Smith Logically Collective on PC 827c2be2410SBarry Smith 828c2be2410SBarry Smith Input Parameters: 8293fc8bf9cSBarry Smith . pc - the multigrid context 830c2be2410SBarry Smith 831c2be2410SBarry Smith Options Database Key: 832c2be2410SBarry Smith $ -pc_mg_galerkin 833c2be2410SBarry Smith 834c2be2410SBarry Smith Level: intermediate 835c2be2410SBarry Smith 836c2be2410SBarry Smith .keywords: MG, set, Galerkin 837c2be2410SBarry Smith 8383fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 8393fc8bf9cSBarry Smith 840c2be2410SBarry Smith @*/ 8417087cfbeSBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc) 842c2be2410SBarry Smith { 843f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 844c2be2410SBarry Smith 845c2be2410SBarry Smith PetscFunctionBegin; 8460700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 847218a07d4SBarry Smith mg->galerkin = PETSC_TRUE; 848c2be2410SBarry Smith PetscFunctionReturn(0); 849c2be2410SBarry Smith } 850c2be2410SBarry Smith 851c2be2410SBarry Smith #undef __FUNCT__ 8523fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 8533fc8bf9cSBarry Smith /*@ 8543fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 8553fc8bf9cSBarry Smith A_i-1 = r_i * A_i * r_i^t 8563fc8bf9cSBarry Smith 8573fc8bf9cSBarry Smith Not Collective 8583fc8bf9cSBarry Smith 8593fc8bf9cSBarry Smith Input Parameter: 8603fc8bf9cSBarry Smith . pc - the multigrid context 8613fc8bf9cSBarry Smith 8623fc8bf9cSBarry Smith Output Parameter: 8633fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 8643fc8bf9cSBarry Smith 8653fc8bf9cSBarry Smith Options Database Key: 8663fc8bf9cSBarry Smith $ -pc_mg_galerkin 8673fc8bf9cSBarry Smith 8683fc8bf9cSBarry Smith Level: intermediate 8693fc8bf9cSBarry Smith 8703fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 8713fc8bf9cSBarry Smith 8723fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 8733fc8bf9cSBarry Smith 8743fc8bf9cSBarry Smith @*/ 8757087cfbeSBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PetscBool *galerkin) 8763fc8bf9cSBarry Smith { 877f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8783fc8bf9cSBarry Smith 8793fc8bf9cSBarry Smith PetscFunctionBegin; 8800700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 881218a07d4SBarry Smith *galerkin = mg->galerkin; 8823fc8bf9cSBarry Smith PetscFunctionReturn(0); 8833fc8bf9cSBarry Smith } 8843fc8bf9cSBarry Smith 8853fc8bf9cSBarry Smith #undef __FUNCT__ 8869dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 8874b9ad928SBarry Smith /*@ 88897177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 88997177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 8904b9ad928SBarry Smith pre-smoothing steps on different levels. 8914b9ad928SBarry Smith 892ad4df100SBarry Smith Logically Collective on PC 8934b9ad928SBarry Smith 8944b9ad928SBarry Smith Input Parameters: 8954b9ad928SBarry Smith + mg - the multigrid context 8964b9ad928SBarry Smith - n - the number of smoothing steps 8974b9ad928SBarry Smith 8984b9ad928SBarry Smith Options Database Key: 8994b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 9004b9ad928SBarry Smith 9014b9ad928SBarry Smith Level: advanced 9024b9ad928SBarry Smith 9034b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 9044b9ad928SBarry Smith 90597177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 9064b9ad928SBarry Smith @*/ 9077087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothDown(PC pc,PetscInt n) 9084b9ad928SBarry Smith { 909f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 910f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 9116849ba73SBarry Smith PetscErrorCode ierr; 91279416396SBarry Smith PetscInt i,levels; 9134b9ad928SBarry Smith 9144b9ad928SBarry Smith PetscFunctionBegin; 9150700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 916e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 917c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 918f3fbd535SBarry Smith levels = mglevels[0]->levels; 9194b9ad928SBarry Smith 920b05257ddSBarry Smith for (i=1; i<levels; i++) { 9214b9ad928SBarry Smith /* make sure smoother up and down are different */ 92297177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 923f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 92431567311SBarry Smith mg->default_smoothd = n; 9254b9ad928SBarry Smith } 9264b9ad928SBarry Smith PetscFunctionReturn(0); 9274b9ad928SBarry Smith } 9284b9ad928SBarry Smith 9294b9ad928SBarry Smith #undef __FUNCT__ 9309dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 9314b9ad928SBarry Smith /*@ 93297177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 93397177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 9344b9ad928SBarry Smith post-smoothing steps on different levels. 9354b9ad928SBarry Smith 936ad4df100SBarry Smith Logically Collective on PC 9374b9ad928SBarry Smith 9384b9ad928SBarry Smith Input Parameters: 9394b9ad928SBarry Smith + mg - the multigrid context 9404b9ad928SBarry Smith - n - the number of smoothing steps 9414b9ad928SBarry Smith 9424b9ad928SBarry Smith Options Database Key: 9434b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 9444b9ad928SBarry Smith 9454b9ad928SBarry Smith Level: advanced 9464b9ad928SBarry Smith 9474b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 948a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 9494b9ad928SBarry Smith 9504b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 9514b9ad928SBarry Smith 95297177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 9534b9ad928SBarry Smith @*/ 9547087cfbeSBarry Smith PetscErrorCode PCMGSetNumberSmoothUp(PC pc,PetscInt n) 9554b9ad928SBarry Smith { 956f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 957f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 9586849ba73SBarry Smith PetscErrorCode ierr; 95979416396SBarry Smith PetscInt i,levels; 9604b9ad928SBarry Smith 9614b9ad928SBarry Smith PetscFunctionBegin; 9620700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 963e7e72b3dSBarry Smith if (!mglevels) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 964c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 965f3fbd535SBarry Smith levels = mglevels[0]->levels; 9664b9ad928SBarry Smith 9674b9ad928SBarry Smith for (i=1; i<levels; i++) { 9684b9ad928SBarry Smith /* make sure smoother up and down are different */ 96997177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 970f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 97131567311SBarry Smith mg->default_smoothu = n; 9724b9ad928SBarry Smith } 9734b9ad928SBarry Smith PetscFunctionReturn(0); 9744b9ad928SBarry Smith } 9754b9ad928SBarry Smith 9764b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 9774b9ad928SBarry Smith 9783b09bd56SBarry Smith /*MC 979ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 9803b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 9813b09bd56SBarry Smith 9823b09bd56SBarry Smith Options Database Keys: 9833b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 9840d353602SBarry Smith . -pc_mg_cycles v or w 98579416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 9863b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 9873b09bd56SBarry Smith . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 9883b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 9893b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 99068eff7e6SBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators 9913b09bd56SBarry Smith - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 992e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 9933b09bd56SBarry Smith 99424c3aa18SBarry Smith 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 9953b09bd56SBarry Smith 9963b09bd56SBarry Smith Level: intermediate 9973b09bd56SBarry Smith 9988f87f92bSBarry Smith Concepts: multigrid/multilevel 9993b09bd56SBarry Smith 100024c3aa18SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC 10010d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 100297177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 100397177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 10040d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 10053b09bd56SBarry Smith M*/ 10063b09bd56SBarry Smith 10074b9ad928SBarry Smith EXTERN_C_BEGIN 10084b9ad928SBarry Smith #undef __FUNCT__ 10094b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 10107087cfbeSBarry Smith PetscErrorCode PCCreate_MG(PC pc) 10114b9ad928SBarry Smith { 1012f3fbd535SBarry Smith PC_MG *mg; 1013f3fbd535SBarry Smith PetscErrorCode ierr; 1014f3fbd535SBarry Smith 10154b9ad928SBarry Smith PetscFunctionBegin; 1016f3fbd535SBarry Smith ierr = PetscNewLog(pc,PC_MG,&mg);CHKERRQ(ierr); 1017f3fbd535SBarry Smith pc->data = (void*)mg; 1018f3fbd535SBarry Smith mg->nlevels = -1; 1019f3fbd535SBarry Smith 10204b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 10214b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 10224b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 10234b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 10244b9ad928SBarry Smith pc->ops->view = PCView_MG; 10254b9ad928SBarry Smith PetscFunctionReturn(0); 10264b9ad928SBarry Smith } 10274b9ad928SBarry Smith EXTERN_C_END 1028