1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines the multigrid preconditioner interface. 44b9ad928SBarry Smith */ 5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 61e25c274SJed Brown #include <petscdm.h> 74b9ad928SBarry Smith 831567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 94b9ad928SBarry Smith { 1031567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1131567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 126849ba73SBarry Smith PetscErrorCode ierr; 1331567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 14a0808db4SHong Zhang PC subpc; 15a0808db4SHong Zhang PCFailedReason pcreason; 164b9ad928SBarry Smith 174b9ad928SBarry Smith PetscFunctionBegin; 1863e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 1931567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 20a0808db4SHong Zhang ierr = KSPGetPC(mglevels->smoothd,&subpc);CHKERRQ(ierr); 21a0808db4SHong Zhang ierr = PCGetSetUpFailedReason(subpc,&pcreason);CHKERRQ(ierr); 22a0808db4SHong Zhang if (pcreason) { 23a0808db4SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 24a0808db4SHong Zhang } 2563e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2631567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 2763e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 2831567311SBarry Smith ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 2963e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 304b9ad928SBarry Smith 314b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 3231567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 334b9ad928SBarry Smith PetscReal rnorm; 3431567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 354b9ad928SBarry Smith if (rnorm <= mg->ttol) { 3670441072SBarry Smith if (rnorm < mg->abstol) { 374d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 3857622a8eSBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr); 394b9ad928SBarry Smith } else { 404d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 4157622a8eSBarry Smith ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than relative tolerance times initial residual norm %g\n",(double)rnorm,(double)mg->ttol);CHKERRQ(ierr); 424b9ad928SBarry Smith } 434b9ad928SBarry Smith PetscFunctionReturn(0); 444b9ad928SBarry Smith } 454b9ad928SBarry Smith } 464b9ad928SBarry Smith 4731567311SBarry Smith mgc = *(mglevelsin - 1); 4863e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 4931567311SBarry Smith ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 5063e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 51efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 524b9ad928SBarry Smith while (cycles--) { 5331567311SBarry Smith ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 544b9ad928SBarry Smith } 5563e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5631567311SBarry Smith ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 5763e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5863e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 5931567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 6063e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 614b9ad928SBarry Smith } 624b9ad928SBarry Smith PetscFunctionReturn(0); 634b9ad928SBarry Smith } 644b9ad928SBarry Smith 65ace3abfcSBarry 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) 664b9ad928SBarry Smith { 67f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 68f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 69dfbe8321SBarry Smith PetscErrorCode ierr; 70f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 714b9ad928SBarry Smith 724b9ad928SBarry Smith PetscFunctionBegin; 73a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 74a762d673SBarry Smith for (i=0; i<levels; i++) { 75a762d673SBarry Smith if (!mglevels[i]->A) { 76a762d673SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 77a762d673SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 78a762d673SBarry Smith } 79a762d673SBarry Smith } 80f3fbd535SBarry Smith mglevels[levels-1]->b = b; 81f3fbd535SBarry Smith mglevels[levels-1]->x = x; 824b9ad928SBarry Smith 8331567311SBarry Smith mg->rtol = rtol; 8431567311SBarry Smith mg->abstol = abstol; 8531567311SBarry Smith mg->dtol = dtol; 864b9ad928SBarry Smith if (rtol) { 874b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 884b9ad928SBarry Smith PetscReal rnorm; 897319c654SBarry Smith if (zeroguess) { 907319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 917319c654SBarry Smith } else { 92f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 934b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 947319c654SBarry Smith } 9531567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 962fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 972fa5cd67SKarl Rupp else mg->ttol = 0.0; 984b9ad928SBarry Smith 994d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 100336babb1SJed Brown stop prematurely due to small residual */ 1014d0a8057SBarry Smith for (i=1; i<levels; i++) { 102f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 103f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 10423067569SBarry Smith /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 10523067569SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 106f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 1074b9ad928SBarry Smith } 1084d0a8057SBarry Smith } 1094d0a8057SBarry Smith 1104d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1114d0a8057SBarry Smith for (i=0; i<its; i++) { 112f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1134d0a8057SBarry Smith if (*reason) break; 1144d0a8057SBarry Smith } 1154d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1164d0a8057SBarry Smith *outits = i; 1174b9ad928SBarry Smith PetscFunctionReturn(0); 1184b9ad928SBarry Smith } 1194b9ad928SBarry Smith 1203aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 1213aeaf226SBarry Smith { 1223aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1233aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1243aeaf226SBarry Smith PetscErrorCode ierr; 1253aeaf226SBarry Smith PetscInt i,n; 1263aeaf226SBarry Smith 1273aeaf226SBarry Smith PetscFunctionBegin; 1283aeaf226SBarry Smith if (mglevels) { 1293aeaf226SBarry Smith n = mglevels[0]->levels; 1303aeaf226SBarry Smith for (i=0; i<n-1; i++) { 1313aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 1323aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 1333aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 1343aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 1353aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 13673250ac0SBarry Smith ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 1373aeaf226SBarry Smith } 1383aeaf226SBarry Smith 1393aeaf226SBarry Smith for (i=0; i<n; i++) { 1403aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 1413aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1423aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 1433aeaf226SBarry Smith } 1443aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 1453aeaf226SBarry Smith } 1463aeaf226SBarry Smith } 1473aeaf226SBarry Smith PetscFunctionReturn(0); 1483aeaf226SBarry Smith } 1493aeaf226SBarry Smith 1504b9ad928SBarry Smith /*@C 15197177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 1524b9ad928SBarry Smith Must be called before any other MG routine. 1534b9ad928SBarry Smith 154ad4df100SBarry Smith Logically Collective on PC 1554b9ad928SBarry Smith 1564b9ad928SBarry Smith Input Parameters: 1574b9ad928SBarry Smith + pc - the preconditioner context 1584b9ad928SBarry Smith . levels - the number of levels 1594b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 1602bf68e3eSBarry Smith on smaller sets of processors. 1614b9ad928SBarry Smith 1624b9ad928SBarry Smith Level: intermediate 1634b9ad928SBarry Smith 1644b9ad928SBarry Smith Notes: 1654b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 1664b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 1674b9ad928SBarry Smith 1684b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 1694b9ad928SBarry Smith 17097177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 1714b9ad928SBarry Smith @*/ 1727087cfbeSBarry Smith PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 1734b9ad928SBarry Smith { 174dfbe8321SBarry Smith PetscErrorCode ierr; 175f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 176ce94432eSBarry Smith MPI_Comm comm; 1773aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 17810eca3edSLisandro Dalcin PCMGType mgtype = mg->am; 17910167fecSBarry Smith PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V; 180f3fbd535SBarry Smith PetscInt i; 181f3fbd535SBarry Smith PetscMPIInt size; 182f3fbd535SBarry Smith const char *prefix; 183f3fbd535SBarry Smith PC ipc; 1843aeaf226SBarry Smith PetscInt n; 1854b9ad928SBarry Smith 1864b9ad928SBarry Smith PetscFunctionBegin; 1870700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 188548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 189548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 1901a97d7b6SStefano Zampini ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1913aeaf226SBarry Smith if (mglevels) { 19210eca3edSLisandro Dalcin mgctype = mglevels[0]->cycles; 1933aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 1943aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 1953aeaf226SBarry Smith n = mglevels[0]->levels; 1963aeaf226SBarry Smith for (i=0; i<n; i++) { 1973aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1983aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 1993aeaf226SBarry Smith } 2003aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 2013aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 2023aeaf226SBarry Smith } 2033aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 2043aeaf226SBarry Smith } 205f3fbd535SBarry Smith 206f3fbd535SBarry Smith mg->nlevels = levels; 207f3fbd535SBarry Smith 208785e854fSJed Brown ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr); 2093bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 210f3fbd535SBarry Smith 211f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 212f3fbd535SBarry Smith 213a9db87a2SMatthew G Knepley mg->stageApply = 0; 214f3fbd535SBarry Smith for (i=0; i<levels; i++) { 215b00a9115SJed Brown ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr); 2162fa5cd67SKarl Rupp 217f3fbd535SBarry Smith mglevels[i]->level = i; 218f3fbd535SBarry Smith mglevels[i]->levels = levels; 21910eca3edSLisandro Dalcin mglevels[i]->cycles = mgctype; 220336babb1SJed Brown mg->default_smoothu = 2; 221336babb1SJed Brown mg->default_smoothd = 2; 22263e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 22363e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 22463e6d426SJed Brown mglevels[i]->eventresidual = 0; 22563e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 226f3fbd535SBarry Smith 227f3fbd535SBarry Smith if (comms) comm = comms[i]; 228f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 229422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr); 2300ee9a668SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 2310ee9a668SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 2320ee9a668SBarry Smith ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 2330ee9a668SBarry Smith if (i || levels == 1) { 2340ee9a668SBarry Smith char tprefix[128]; 2350ee9a668SBarry Smith 236336babb1SJed Brown ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 2370059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 238336babb1SJed Brown ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 239336babb1SJed Brown ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 240336babb1SJed Brown ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 2410ee9a668SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 242f3fbd535SBarry Smith 2430ee9a668SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 2440ee9a668SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 2450ee9a668SBarry Smith } else { 246f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 247f3fbd535SBarry Smith 2489dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 249f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 250f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 251f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 252f3fbd535SBarry Smith if (size > 1) { 253f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 254f3fbd535SBarry Smith } else { 255f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 256f3fbd535SBarry Smith } 257753b7fb9SBarry Smith ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 258f3fbd535SBarry Smith } 2593bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr); 2602fa5cd67SKarl Rupp 261f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 26231567311SBarry Smith mg->rtol = 0.0; 26331567311SBarry Smith mg->abstol = 0.0; 26431567311SBarry Smith mg->dtol = 0.0; 26531567311SBarry Smith mg->ttol = 0.0; 26631567311SBarry Smith mg->cyclesperpcapply = 1; 267f3fbd535SBarry Smith } 268f3fbd535SBarry Smith mg->levels = mglevels; 26910eca3edSLisandro Dalcin ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 2704b9ad928SBarry Smith PetscFunctionReturn(0); 2714b9ad928SBarry Smith } 2724b9ad928SBarry Smith 273c07bf074SBarry Smith 274c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 275c07bf074SBarry Smith { 276c07bf074SBarry Smith PetscErrorCode ierr; 277a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 278a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 279a06653b4SBarry Smith PetscInt i,n; 280c07bf074SBarry Smith 281c07bf074SBarry Smith PetscFunctionBegin; 282a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 283a06653b4SBarry Smith if (mglevels) { 284a06653b4SBarry Smith n = mglevels[0]->levels; 285a06653b4SBarry Smith for (i=0; i<n; i++) { 286a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2876bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 288a06653b4SBarry Smith } 2896bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 290a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 291a06653b4SBarry Smith } 292a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 293a06653b4SBarry Smith } 294c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 295f3fbd535SBarry Smith PetscFunctionReturn(0); 296f3fbd535SBarry Smith } 297f3fbd535SBarry Smith 298f3fbd535SBarry Smith 299f3fbd535SBarry Smith 30009573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 30109573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 30209573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 303f3fbd535SBarry Smith 304f3fbd535SBarry Smith /* 305f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 306f3fbd535SBarry Smith or full cycle of multigrid. 307f3fbd535SBarry Smith 308f3fbd535SBarry Smith Note: 309f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 310f3fbd535SBarry Smith */ 311f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 312f3fbd535SBarry Smith { 313f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 314f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 315f3fbd535SBarry Smith PetscErrorCode ierr; 316f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 317f3fbd535SBarry Smith 318f3fbd535SBarry Smith PetscFunctionBegin; 319a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 320e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 321298cc208SBarry Smith for (i=0; i<levels; i++) { 322e1d8e5deSBarry Smith if (!mglevels[i]->A) { 32323ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 324298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 325e1d8e5deSBarry Smith } 326e1d8e5deSBarry Smith } 327e1d8e5deSBarry Smith 328f3fbd535SBarry Smith mglevels[levels-1]->b = b; 329f3fbd535SBarry Smith mglevels[levels-1]->x = x; 33031567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 331f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 33231567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 3330298fd71SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr); 334f3fbd535SBarry Smith } 3352fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 33631567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 3372fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 33831567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 3392fa5cd67SKarl Rupp } else { 340f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 341f3fbd535SBarry Smith } 342a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 343f3fbd535SBarry Smith PetscFunctionReturn(0); 344f3fbd535SBarry Smith } 345f3fbd535SBarry Smith 346f3fbd535SBarry Smith 3474416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 348f3fbd535SBarry Smith { 349f3fbd535SBarry Smith PetscErrorCode ierr; 350710315b6SLawrence Mitchell PetscInt levels,cycles; 3512134b1e4SBarry Smith PetscBool flg; 352f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 3533d3eaba7SBarry Smith PC_MG_Levels **mglevels; 354f3fbd535SBarry Smith PCMGType mgtype; 355f3fbd535SBarry Smith PCMGCycleType mgctype; 3562134b1e4SBarry Smith PCMGGalerkinType gtype; 357f3fbd535SBarry Smith 358f3fbd535SBarry Smith PetscFunctionBegin; 3591d743356SStefano Zampini levels = PetscMax(mg->nlevels,1); 360e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr); 361f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 3621a97d7b6SStefano Zampini if (!flg && !mg->levels && pc->dm) { 363eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 364eb3f98d2SBarry Smith levels++; 3653aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 366eb3f98d2SBarry Smith } 3670298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 3683aeaf226SBarry Smith mglevels = mg->levels; 3693aeaf226SBarry Smith 370f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 371f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 372f3fbd535SBarry Smith if (flg) { 373f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 3742fa5cd67SKarl Rupp } 3752134b1e4SBarry Smith gtype = mg->galerkin; 3762134b1e4SBarry Smith ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg);CHKERRQ(ierr); 3772134b1e4SBarry Smith if (flg) { 3782134b1e4SBarry Smith ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr); 379f3fbd535SBarry Smith } 380f56b1265SBarry Smith flg = PETSC_FALSE; 381f442ab6aSBarry Smith ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create seperate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr); 382f442ab6aSBarry Smith if (flg) { 383f442ab6aSBarry Smith ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr); 384f442ab6aSBarry Smith } 38531567311SBarry Smith mgtype = mg->am; 386f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 387f3fbd535SBarry Smith if (flg) { 388f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 389f3fbd535SBarry Smith } 39031567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 3915363c1e0SLisandro Dalcin ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 392f3fbd535SBarry Smith if (flg) { 393f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 394f3fbd535SBarry Smith } 395f3fbd535SBarry Smith } 396f3fbd535SBarry Smith flg = PETSC_FALSE; 3970298fd71SBarry Smith ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 398f3fbd535SBarry Smith if (flg) { 399f3fbd535SBarry Smith PetscInt i; 400f3fbd535SBarry Smith char eventname[128]; 4011a97d7b6SStefano Zampini 402f3fbd535SBarry Smith levels = mglevels[0]->levels; 403f3fbd535SBarry Smith for (i=0; i<levels; i++) { 404f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 40563e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 406f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 40763e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 408f3fbd535SBarry Smith if (i) { 409f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 41063e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 411f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 41263e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 413f3fbd535SBarry Smith } 414f3fbd535SBarry Smith } 415eec35531SMatthew G Knepley 416ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 417eec35531SMatthew G Knepley { 418eec35531SMatthew G Knepley const char *sname = "MG Apply"; 419eec35531SMatthew G Knepley PetscStageLog stageLog; 420eec35531SMatthew G Knepley PetscInt st; 421eec35531SMatthew G Knepley 422eec35531SMatthew G Knepley PetscFunctionBegin; 423eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 424eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 425eec35531SMatthew G Knepley PetscBool same; 426eec35531SMatthew G Knepley 427eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 4282fa5cd67SKarl Rupp if (same) mg->stageApply = st; 429eec35531SMatthew G Knepley } 430eec35531SMatthew G Knepley if (!mg->stageApply) { 431eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 432eec35531SMatthew G Knepley } 433eec35531SMatthew G Knepley } 434ec60ca73SSatish Balay #endif 435f3fbd535SBarry Smith } 436f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 437f3fbd535SBarry Smith PetscFunctionReturn(0); 438f3fbd535SBarry Smith } 439f3fbd535SBarry Smith 4406a6fc655SJed Brown const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 4416a6fc655SJed Brown const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 4422134b1e4SBarry Smith const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",0}; 443f3fbd535SBarry Smith 4449804daf3SBarry Smith #include <petscdraw.h> 445f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 446f3fbd535SBarry Smith { 447f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 448f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 449f3fbd535SBarry Smith PetscErrorCode ierr; 450e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 451219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 452f3fbd535SBarry Smith 453f3fbd535SBarry Smith PetscFunctionBegin; 454251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 4555b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 456219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 457f3fbd535SBarry Smith if (iascii) { 458e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 459efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 46031567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 46131567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 462f3fbd535SBarry Smith } 4632134b1e4SBarry Smith if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 464f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 4652134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 4662134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr); 4672134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 4682134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr); 4692134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 4702134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr); 4714f66f45eSBarry Smith } else { 4724f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 473f3fbd535SBarry Smith } 4745adeb434SBarry Smith if (mg->view){ 4755adeb434SBarry Smith ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr); 4765adeb434SBarry Smith } 477f3fbd535SBarry Smith for (i=0; i<levels; i++) { 478f3fbd535SBarry Smith if (!i) { 47989cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 480f3fbd535SBarry Smith } else { 48189cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 482f3fbd535SBarry Smith } 483f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 484f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 485f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 486f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 487f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 488f3fbd535SBarry Smith } else if (i) { 48989cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 490f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 491f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 492f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 493f3fbd535SBarry Smith } 494f3fbd535SBarry Smith } 4955b0b0462SBarry Smith } else if (isbinary) { 4965b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 4975b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 4985b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 4995b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 5005b0b0462SBarry Smith } 5015b0b0462SBarry Smith } 502219b1045SBarry Smith } else if (isdraw) { 503219b1045SBarry Smith PetscDraw draw; 504b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 505219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 506219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 5070298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 508219b1045SBarry Smith bottom = y - th; 509219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 510b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 511219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 512219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 513219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 514b4375e8dSPeter Brune } else { 515b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 516b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 517b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 518b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 519b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 520b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 521b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 522b4375e8dSPeter Brune } 5230298fd71SBarry Smith ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 5241cd381d6SBarry Smith bottom -= th; 525219b1045SBarry Smith } 5265b0b0462SBarry Smith } 527f3fbd535SBarry Smith PetscFunctionReturn(0); 528f3fbd535SBarry Smith } 529f3fbd535SBarry Smith 530af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 531af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 532cab2e9ccSBarry Smith 533f3fbd535SBarry Smith /* 534f3fbd535SBarry Smith Calls setup for the KSP on each level 535f3fbd535SBarry Smith */ 536f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 537f3fbd535SBarry Smith { 538f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 539f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 540f3fbd535SBarry Smith PetscErrorCode ierr; 5411a97d7b6SStefano Zampini PetscInt i,n; 54298e04984SBarry Smith PC cpc; 5433db492dfSBarry Smith PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 544f3fbd535SBarry Smith Mat dA,dB; 545f3fbd535SBarry Smith Vec tvec; 546218a07d4SBarry Smith DM *dms; 547649052a6SBarry Smith PetscViewer viewer = 0; 5482134b1e4SBarry Smith PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE; 549f3fbd535SBarry Smith 550f3fbd535SBarry Smith PetscFunctionBegin; 5511a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up"); 5521a97d7b6SStefano Zampini n = mglevels[0]->levels; 55301bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 5543aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 5553aeaf226SBarry Smith PetscInt levels; 5563aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 5573aeaf226SBarry Smith levels++; 5583aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 5590298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 5603aeaf226SBarry Smith n = levels; 5613aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 5623aeaf226SBarry Smith mglevels = mg->levels; 5633aeaf226SBarry Smith } 5643aeaf226SBarry Smith } 56598e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 5663aeaf226SBarry Smith 567f3fbd535SBarry Smith 568f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 569f3fbd535SBarry Smith /* so use those from global PC */ 570f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 5710298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 57298e04984SBarry Smith if (opsset) { 57398e04984SBarry Smith Mat mmat; 57423ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 57598e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 57698e04984SBarry Smith } 5774a5f07a5SMark F. Adams 5784a5f07a5SMark F. Adams if (!opsset) { 57971b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 5804a5f07a5SMark F. Adams if(use_amat){ 581f3fbd535SBarry 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); 58223ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 583f3fbd535SBarry Smith } 5844a5f07a5SMark F. Adams else { 5854a5f07a5SMark F. Adams ierr = PetscInfo(pc,"Using matrix (pmat) operators to define finest grid operator \n because PCMGGetSmoother(pc,nlevels-1,&ksp);KSPSetOperators(ksp,...); was not called.\n");CHKERRQ(ierr); 58623ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 5874a5f07a5SMark F. Adams } 5884a5f07a5SMark F. Adams } 589f3fbd535SBarry Smith 5909c7ffb39SBarry Smith for (i=n-1; i>0; i--) { 5919c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 5929c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 5939c7ffb39SBarry Smith continue; 5949c7ffb39SBarry Smith } 5959c7ffb39SBarry Smith } 5962134b1e4SBarry Smith 5972134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 5982134b1e4SBarry Smith if (dA == dB) dAeqdB = PETSC_TRUE; 5992134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) { 6002134b1e4SBarry Smith needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 6012134b1e4SBarry Smith } 6022134b1e4SBarry Smith 6032134b1e4SBarry Smith 6049c7ffb39SBarry Smith /* 6059c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 6069c7ffb39SBarry Smith Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 6079c7ffb39SBarry Smith */ 6082134b1e4SBarry Smith if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 6092d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 6102e499ae9SBarry Smith Mat p; 61173250ac0SBarry Smith Vec rscale; 612785e854fSJed Brown ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 613218a07d4SBarry Smith dms[n-1] = pc->dm; 614ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 615ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 61641fce8fdSHong Zhang /* 61741fce8fdSHong Zhang Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level. 61841fce8fdSHong Zhang Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options. 61941fce8fdSHong Zhang But it is safe to use -dm_mat_type. 62041fce8fdSHong Zhang 62141fce8fdSHong Zhang The mat type should not be hardcoded like this, we need to find a better way. 62241fce8fdSHong Zhang ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr); 62341fce8fdSHong Zhang */ 624218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 625942e3340SBarry Smith DMKSP kdm; 626*eab52d2bSLawrence Mitchell PetscBool dmhasrestrict, dmhasinject; 6273c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 6282134b1e4SBarry Smith if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 629942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 630d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 631d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 6320298fd71SBarry Smith kdm->ops->computerhs = NULL; 6330298fd71SBarry Smith kdm->rhsctx = NULL; 63424c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 635e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 6362d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 63700ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 63873250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 6396bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 640218a07d4SBarry Smith } 6413ad4599aSBarry Smith ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr); 6423ad4599aSBarry Smith if (dmhasrestrict && !mglevels[i+1]->restrct){ 6433ad4599aSBarry Smith ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr); 6443ad4599aSBarry Smith ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr); 6453ad4599aSBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 6463ad4599aSBarry Smith } 647*eab52d2bSLawrence Mitchell ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr); 648*eab52d2bSLawrence Mitchell if (dmhasinject && !mglevels[i+1]->inject){ 649*eab52d2bSLawrence Mitchell ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr); 650*eab52d2bSLawrence Mitchell ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr); 651*eab52d2bSLawrence Mitchell ierr = MatDestroy(&p);CHKERRQ(ierr); 652*eab52d2bSLawrence Mitchell } 65324c3aa18SBarry Smith } 6542d2b81a6SBarry Smith 655ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 6562d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 657ce4cda84SJed Brown } 658cab2e9ccSBarry Smith 659ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 6602134b1e4SBarry Smith /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 661cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 662cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 663218a07d4SBarry Smith } 664218a07d4SBarry Smith 6652134b1e4SBarry Smith if (mg->galerkin < PC_MG_GALERKIN_NONE) { 6662134b1e4SBarry Smith Mat A,B; 6672134b1e4SBarry Smith PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 6682134b1e4SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 6692134b1e4SBarry Smith 6702134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 6712134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE; 6722134b1e4SBarry Smith if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 673f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 674b9367914SBarry Smith if (!mglevels[i+1]->restrct && !mglevels[i+1]->interpolate) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must provide interpolation or restriction for each MG level except level 0"); 675b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 676b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 677b9367914SBarry Smith } 678b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 679b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 680b9367914SBarry Smith } 6812134b1e4SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 6822134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr); 6832134b1e4SBarry Smith } 6842134b1e4SBarry Smith if (doA) { 6852df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr); 6862134b1e4SBarry Smith } 6872134b1e4SBarry Smith if (doB) { 6882df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr); 689b9367914SBarry Smith } 6902134b1e4SBarry Smith /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 6912134b1e4SBarry Smith if (!doA && dAeqdB) { 6922134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);} 6932134b1e4SBarry Smith A = B; 6942134b1e4SBarry Smith } else if (!doA && reuse == MAT_INITIAL_MATRIX ) { 6952134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr); 6962134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 697b9367914SBarry Smith } 6982134b1e4SBarry Smith if (!doB && dAeqdB) { 6992134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);} 7002134b1e4SBarry Smith B = A; 7012134b1e4SBarry Smith } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 70223ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 7032134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 7047d537d92SJed Brown } 7052134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 7062134b1e4SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr); 7072134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr); 7082134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr); 7092134b1e4SBarry Smith } 7102134b1e4SBarry Smith dA = A; 711f3fbd535SBarry Smith dB = B; 712f3fbd535SBarry Smith } 713f3fbd535SBarry Smith } 7142134b1e4SBarry Smith if (needRestricts && pc->dm && pc->dm->x) { 715cab2e9ccSBarry Smith /* need to restrict Jacobian location to coarser meshes for evaluation */ 716cab2e9ccSBarry Smith for (i=n-2; i>-1; i--) { 717c88c5224SJed Brown Mat R; 718c88c5224SJed Brown Vec rscale; 719cab2e9ccSBarry Smith if (!mglevels[i]->smoothd->dm->x) { 720cab2e9ccSBarry Smith Vec *vecs; 7212a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vecs,0,NULL);CHKERRQ(ierr); 722cab2e9ccSBarry Smith mglevels[i]->smoothd->dm->x = vecs[0]; 723cab2e9ccSBarry Smith ierr = PetscFree(vecs);CHKERRQ(ierr); 724cab2e9ccSBarry Smith } 725c88c5224SJed Brown ierr = PCMGGetRestriction(pc,i+1,&R);CHKERRQ(ierr); 726c88c5224SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 727c88c5224SJed Brown ierr = MatRestrict(R,mglevels[i+1]->smoothd->dm->x,mglevels[i]->smoothd->dm->x);CHKERRQ(ierr); 728c88c5224SJed Brown ierr = VecPointwiseMult(mglevels[i]->smoothd->dm->x,mglevels[i]->smoothd->dm->x,rscale);CHKERRQ(ierr); 729cab2e9ccSBarry Smith } 730f3fbd535SBarry Smith } 7312134b1e4SBarry Smith if (needRestricts && pc->dm) { 732caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 733ccceb50cSJed Brown DM dmfine,dmcoarse; 734caa4e7f2SJed Brown Mat Restrict,Inject; 735caa4e7f2SJed Brown Vec rscale; 736ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 737ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 738caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 739caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 740*eab52d2bSLawrence Mitchell ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr); 741caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 742caa4e7f2SJed Brown } 743caa4e7f2SJed Brown } 744f3fbd535SBarry Smith 745f3fbd535SBarry Smith if (!pc->setupcalled) { 746f3fbd535SBarry Smith for (i=0; i<n; i++) { 747f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 748f3fbd535SBarry Smith } 749f3fbd535SBarry Smith for (i=1; i<n; i++) { 750f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 751f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 752f3fbd535SBarry Smith } 753f3fbd535SBarry Smith } 7543ad4599aSBarry Smith /* insure that if either interpolation or restriction is set the other other one is set */ 755f3fbd535SBarry Smith for (i=1; i<n; i++) { 7563ad4599aSBarry Smith ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr); 7573ad4599aSBarry Smith ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr); 758f3fbd535SBarry Smith } 759f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 760f3fbd535SBarry Smith if (!mglevels[i]->b) { 761f3fbd535SBarry Smith Vec *vec; 7622a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 763f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 7646bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 765f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 766f3fbd535SBarry Smith } 767f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 768f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 769f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 7706bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 771f3fbd535SBarry Smith } 772f3fbd535SBarry Smith if (!mglevels[i]->x) { 773f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 774f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 7756bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 776f3fbd535SBarry Smith } 777f3fbd535SBarry Smith } 778f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 779f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 780f3fbd535SBarry Smith Vec *vec; 7812a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 782f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 7836bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 784f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 785f3fbd535SBarry Smith } 786f3fbd535SBarry Smith } 787f3fbd535SBarry Smith 78898e04984SBarry Smith if (pc->dm) { 78998e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 79098e04984SBarry Smith for (i=0; i<n-1; i++) { 79198e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 79298e04984SBarry Smith } 79398e04984SBarry Smith } 794f3fbd535SBarry Smith 795f3fbd535SBarry Smith for (i=1; i<n; i++) { 7962a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 797f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 798f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 799f3fbd535SBarry Smith } 80063e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 801f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 802899639b0SHong Zhang if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 803899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 804899639b0SHong Zhang } 80563e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 806d42688cbSBarry Smith if (!mglevels[i]->residual) { 807d42688cbSBarry Smith Mat mat; 80813842ffbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 80954b2cd4bSJed Brown ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 810d42688cbSBarry Smith } 811f3fbd535SBarry Smith } 812f3fbd535SBarry Smith for (i=1; i<n; i++) { 813f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 814f3fbd535SBarry Smith Mat downmat,downpmat; 815f3fbd535SBarry Smith 816f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 8170298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 818f3fbd535SBarry Smith if (!opsset) { 81923ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 82023ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 821f3fbd535SBarry Smith } 822f3fbd535SBarry Smith 823f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 82463e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 825f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 826899639b0SHong Zhang if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PCSETUP_FAILED) { 827899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 828899639b0SHong Zhang } 82963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 830f3fbd535SBarry Smith } 831f3fbd535SBarry Smith } 832f3fbd535SBarry Smith 83363e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 834f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 835899639b0SHong Zhang if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PCSETUP_FAILED) { 836899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 837899639b0SHong Zhang } 83863e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 839f3fbd535SBarry Smith 840f3fbd535SBarry Smith /* 841f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 842e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 843f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 844f3fbd535SBarry Smith 845f3fbd535SBarry Smith Only support one or the other at the same time. 846f3fbd535SBarry Smith */ 847f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 848c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 849ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 850f3fbd535SBarry Smith dump = PETSC_FALSE; 851f3fbd535SBarry Smith #endif 852c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 853ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 854f3fbd535SBarry Smith 855f3fbd535SBarry Smith if (viewer) { 856f3fbd535SBarry Smith for (i=1; i<n; i++) { 857f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 858f3fbd535SBarry Smith } 859f3fbd535SBarry Smith for (i=0; i<n; i++) { 860f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 861f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 862f3fbd535SBarry Smith } 863f3fbd535SBarry Smith } 864f3fbd535SBarry Smith PetscFunctionReturn(0); 865f3fbd535SBarry Smith } 866f3fbd535SBarry Smith 867f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 868f3fbd535SBarry Smith 8694b9ad928SBarry Smith /*@ 87097177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 8714b9ad928SBarry Smith 8724b9ad928SBarry Smith Not Collective 8734b9ad928SBarry Smith 8744b9ad928SBarry Smith Input Parameter: 8754b9ad928SBarry Smith . pc - the preconditioner context 8764b9ad928SBarry Smith 8774b9ad928SBarry Smith Output parameter: 8784b9ad928SBarry Smith . levels - the number of levels 8794b9ad928SBarry Smith 8804b9ad928SBarry Smith Level: advanced 8814b9ad928SBarry Smith 8824b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 8834b9ad928SBarry Smith 88497177400SBarry Smith .seealso: PCMGSetLevels() 8854b9ad928SBarry Smith @*/ 8867087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 8874b9ad928SBarry Smith { 888f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 8894b9ad928SBarry Smith 8904b9ad928SBarry Smith PetscFunctionBegin; 8910700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8924482741eSBarry Smith PetscValidIntPointer(levels,2); 893f3fbd535SBarry Smith *levels = mg->nlevels; 8944b9ad928SBarry Smith PetscFunctionReturn(0); 8954b9ad928SBarry Smith } 8964b9ad928SBarry Smith 8974b9ad928SBarry Smith /*@ 89897177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 8994b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 9004b9ad928SBarry Smith 901ad4df100SBarry Smith Logically Collective on PC 9024b9ad928SBarry Smith 9034b9ad928SBarry Smith Input Parameters: 9044b9ad928SBarry Smith + pc - the preconditioner context 9059dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 9069dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 9074b9ad928SBarry Smith 9084b9ad928SBarry Smith Options Database Key: 9094b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 9104b9ad928SBarry Smith additive, full, kaskade 9114b9ad928SBarry Smith 9124b9ad928SBarry Smith Level: advanced 9134b9ad928SBarry Smith 9144b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 9154b9ad928SBarry Smith 91697177400SBarry Smith .seealso: PCMGSetLevels() 9174b9ad928SBarry Smith @*/ 9187087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 9194b9ad928SBarry Smith { 920f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 9214b9ad928SBarry Smith 9224b9ad928SBarry Smith PetscFunctionBegin; 9230700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 924c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 92531567311SBarry Smith mg->am = form; 9269dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 927c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 928c60c7ad4SBarry Smith PetscFunctionReturn(0); 929c60c7ad4SBarry Smith } 930c60c7ad4SBarry Smith 931c60c7ad4SBarry Smith /*@ 932c60c7ad4SBarry Smith PCMGGetType - Determines the form of multigrid to use: 933c60c7ad4SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 934c60c7ad4SBarry Smith 935c60c7ad4SBarry Smith Logically Collective on PC 936c60c7ad4SBarry Smith 937c60c7ad4SBarry Smith Input Parameter: 938c60c7ad4SBarry Smith . pc - the preconditioner context 939c60c7ad4SBarry Smith 940c60c7ad4SBarry Smith Output Parameter: 941c60c7ad4SBarry Smith . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 942c60c7ad4SBarry Smith 943c60c7ad4SBarry Smith 944c60c7ad4SBarry Smith Level: advanced 945c60c7ad4SBarry Smith 946c60c7ad4SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 947c60c7ad4SBarry Smith 948c60c7ad4SBarry Smith .seealso: PCMGSetLevels() 949c60c7ad4SBarry Smith @*/ 950c60c7ad4SBarry Smith PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 951c60c7ad4SBarry Smith { 952c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 953c60c7ad4SBarry Smith 954c60c7ad4SBarry Smith PetscFunctionBegin; 955c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 956c60c7ad4SBarry Smith *type = mg->am; 9574b9ad928SBarry Smith PetscFunctionReturn(0); 9584b9ad928SBarry Smith } 9594b9ad928SBarry Smith 9604b9ad928SBarry Smith /*@ 9610d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 9624b9ad928SBarry Smith complicated cycling. 9634b9ad928SBarry Smith 964ad4df100SBarry Smith Logically Collective on PC 9654b9ad928SBarry Smith 9664b9ad928SBarry Smith Input Parameters: 967c2be2410SBarry Smith + pc - the multigrid context 968c1cbb1deSBarry Smith - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 9694b9ad928SBarry Smith 9704b9ad928SBarry Smith Options Database Key: 971c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired 9724b9ad928SBarry Smith 9734b9ad928SBarry Smith Level: advanced 9744b9ad928SBarry Smith 9754b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 9764b9ad928SBarry Smith 9770d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 9784b9ad928SBarry Smith @*/ 9797087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 9804b9ad928SBarry Smith { 981f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 982f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 98379416396SBarry Smith PetscInt i,levels; 9844b9ad928SBarry Smith 9854b9ad928SBarry Smith PetscFunctionBegin; 9860700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 98769fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc,n,2); 9881a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 989f3fbd535SBarry Smith levels = mglevels[0]->levels; 9902fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 9914b9ad928SBarry Smith PetscFunctionReturn(0); 9924b9ad928SBarry Smith } 9934b9ad928SBarry Smith 9948cc2d5dfSBarry Smith /*@ 9958cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 9968cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 9978cc2d5dfSBarry Smith 998ad4df100SBarry Smith Logically Collective on PC 9998cc2d5dfSBarry Smith 10008cc2d5dfSBarry Smith Input Parameters: 10018cc2d5dfSBarry Smith + pc - the multigrid context 10028cc2d5dfSBarry Smith - n - number of cycles (default is 1) 10038cc2d5dfSBarry Smith 10048cc2d5dfSBarry Smith Options Database Key: 1005e1bc860dSBarry Smith . -pc_mg_multiplicative_cycles n 10068cc2d5dfSBarry Smith 10078cc2d5dfSBarry Smith Level: advanced 10088cc2d5dfSBarry Smith 10098cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 10108cc2d5dfSBarry Smith 10118cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 10128cc2d5dfSBarry Smith 10138cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 10148cc2d5dfSBarry Smith @*/ 10157087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 10168cc2d5dfSBarry Smith { 1017f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10188cc2d5dfSBarry Smith 10198cc2d5dfSBarry Smith PetscFunctionBegin; 10200700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1021c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 10222a21e185SBarry Smith mg->cyclesperpcapply = n; 10238cc2d5dfSBarry Smith PetscFunctionReturn(0); 10248cc2d5dfSBarry Smith } 10258cc2d5dfSBarry Smith 10262134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1027fb15c04eSBarry Smith { 1028fb15c04eSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1029fb15c04eSBarry Smith 1030fb15c04eSBarry Smith PetscFunctionBegin; 10312134b1e4SBarry Smith mg->galerkin = use; 1032fb15c04eSBarry Smith PetscFunctionReturn(0); 1033fb15c04eSBarry Smith } 1034fb15c04eSBarry Smith 1035c2be2410SBarry Smith /*@ 103697177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 103782b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1038c2be2410SBarry Smith 1039ad4df100SBarry Smith Logically Collective on PC 1040c2be2410SBarry Smith 1041c2be2410SBarry Smith Input Parameters: 1042c91913d3SJed Brown + pc - the multigrid context 10432134b1e4SBarry Smith - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1044c2be2410SBarry Smith 1045c2be2410SBarry Smith Options Database Key: 10462134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> 1047c2be2410SBarry Smith 1048c2be2410SBarry Smith Level: intermediate 1049c2be2410SBarry Smith 10501c1aac46SBarry Smith Notes: Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 10511c1aac46SBarry Smith use the PCMG construction of the coarser grids. 10521c1aac46SBarry Smith 1053c2be2410SBarry Smith .keywords: MG, set, Galerkin 1054c2be2410SBarry Smith 10552134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType 10563fc8bf9cSBarry Smith 1057c2be2410SBarry Smith @*/ 10582134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1059c2be2410SBarry Smith { 1060fb15c04eSBarry Smith PetscErrorCode ierr; 1061c2be2410SBarry Smith 1062c2be2410SBarry Smith PetscFunctionBegin; 10630700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10642134b1e4SBarry Smith ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr); 1065c2be2410SBarry Smith PetscFunctionReturn(0); 1066c2be2410SBarry Smith } 1067c2be2410SBarry Smith 10683fc8bf9cSBarry Smith /*@ 10693fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 107082b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 10713fc8bf9cSBarry Smith 10723fc8bf9cSBarry Smith Not Collective 10733fc8bf9cSBarry Smith 10743fc8bf9cSBarry Smith Input Parameter: 10753fc8bf9cSBarry Smith . pc - the multigrid context 10763fc8bf9cSBarry Smith 10773fc8bf9cSBarry Smith Output Parameter: 10782134b1e4SBarry Smith . galerkin - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, or PC_MG_GALERKIN_EXTERNAL 10793fc8bf9cSBarry Smith 10803fc8bf9cSBarry Smith Level: intermediate 10813fc8bf9cSBarry Smith 10823fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 10833fc8bf9cSBarry Smith 10842134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType 10853fc8bf9cSBarry Smith 10863fc8bf9cSBarry Smith @*/ 10872134b1e4SBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 10883fc8bf9cSBarry Smith { 1089f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10903fc8bf9cSBarry Smith 10913fc8bf9cSBarry Smith PetscFunctionBegin; 10920700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10932134b1e4SBarry Smith *galerkin = mg->galerkin; 10943fc8bf9cSBarry Smith PetscFunctionReturn(0); 10953fc8bf9cSBarry Smith } 10963fc8bf9cSBarry Smith 10974b9ad928SBarry Smith /*@ 109806792cafSBarry Smith PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1099710315b6SLawrence Mitchell on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1100710315b6SLawrence Mitchell pre- and post-smoothing steps. 110106792cafSBarry Smith 110206792cafSBarry Smith Logically Collective on PC 110306792cafSBarry Smith 110406792cafSBarry Smith Input Parameters: 110506792cafSBarry Smith + mg - the multigrid context 110606792cafSBarry Smith - n - the number of smoothing steps 110706792cafSBarry Smith 110806792cafSBarry Smith Options Database Key: 110906792cafSBarry Smith + -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 111006792cafSBarry Smith 111106792cafSBarry Smith Level: advanced 111206792cafSBarry Smith 111306792cafSBarry Smith Notes: this does not set a value on the coarsest grid, since we assume that 111406792cafSBarry Smith there is no separate smooth up on the coarsest grid. 111506792cafSBarry Smith 111606792cafSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 111706792cafSBarry Smith 1118710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp() 111906792cafSBarry Smith @*/ 112006792cafSBarry Smith PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 112106792cafSBarry Smith { 112206792cafSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 112306792cafSBarry Smith PC_MG_Levels **mglevels = mg->levels; 112406792cafSBarry Smith PetscErrorCode ierr; 112506792cafSBarry Smith PetscInt i,levels; 112606792cafSBarry Smith 112706792cafSBarry Smith PetscFunctionBegin; 112806792cafSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 112906792cafSBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 11301a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 113106792cafSBarry Smith levels = mglevels[0]->levels; 113206792cafSBarry Smith 113306792cafSBarry Smith for (i=1; i<levels; i++) { 113406792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 113506792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 113606792cafSBarry Smith mg->default_smoothu = n; 113706792cafSBarry Smith mg->default_smoothd = n; 113806792cafSBarry Smith } 113906792cafSBarry Smith PetscFunctionReturn(0); 114006792cafSBarry Smith } 114106792cafSBarry Smith 1142f442ab6aSBarry Smith /*@ 1143f442ab6aSBarry Smith PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a seperate KSP from the down (pre) smoother on all levels 1144710315b6SLawrence Mitchell and adds the suffix _up to the options name 1145f442ab6aSBarry Smith 1146f442ab6aSBarry Smith Logically Collective on PC 1147f442ab6aSBarry Smith 1148f442ab6aSBarry Smith Input Parameters: 1149f442ab6aSBarry Smith . pc - the preconditioner context 1150f442ab6aSBarry Smith 1151f442ab6aSBarry Smith Options Database Key: 1152f442ab6aSBarry Smith . -pc_mg_distinct_smoothup 1153f442ab6aSBarry Smith 1154f442ab6aSBarry Smith Level: advanced 1155f442ab6aSBarry Smith 1156f442ab6aSBarry Smith Notes: this does not set a value on the coarsest grid, since we assume that 1157f442ab6aSBarry Smith there is no separate smooth up on the coarsest grid. 1158f442ab6aSBarry Smith 1159f442ab6aSBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 1160f442ab6aSBarry Smith 1161710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth() 1162f442ab6aSBarry Smith @*/ 1163f442ab6aSBarry Smith PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1164f442ab6aSBarry Smith { 1165f442ab6aSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1166f442ab6aSBarry Smith PC_MG_Levels **mglevels = mg->levels; 1167f442ab6aSBarry Smith PetscErrorCode ierr; 1168f442ab6aSBarry Smith PetscInt i,levels; 1169f442ab6aSBarry Smith KSP subksp; 1170f442ab6aSBarry Smith 1171f442ab6aSBarry Smith PetscFunctionBegin; 1172f442ab6aSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1173f442ab6aSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1174f442ab6aSBarry Smith levels = mglevels[0]->levels; 1175f442ab6aSBarry Smith 1176f442ab6aSBarry Smith for (i=1; i<levels; i++) { 1177710315b6SLawrence Mitchell const char *prefix = NULL; 1178f442ab6aSBarry Smith /* make sure smoother up and down are different */ 1179f442ab6aSBarry Smith ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr); 1180710315b6SLawrence Mitchell ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr); 1181710315b6SLawrence Mitchell ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr); 1182f442ab6aSBarry Smith ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr); 1183f442ab6aSBarry Smith } 1184f442ab6aSBarry Smith PetscFunctionReturn(0); 1185f442ab6aSBarry Smith } 1186f442ab6aSBarry Smith 11874b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 11884b9ad928SBarry Smith 11893b09bd56SBarry Smith /*MC 1190ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 11913b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 11923b09bd56SBarry Smith 11933b09bd56SBarry Smith Options Database Keys: 11943b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 11954afba7b0SPatrick Sanan . -pc_mg_cycle_type <v,w> - 11968c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 11973b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 1198710315b6SLawrence Mitchell . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 11992134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 12008cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 12018cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1202e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 12038cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 12048cf6037eSBarry Smith to the binary output file called binaryoutput 12053b09bd56SBarry Smith 1206e941fc33SBarry Smith Notes: If one uses a Krylov method such GMRES or CG as the smoother than one must use KSPFGMRES, KSPGCG, or KSPRICHARDSON as the outer Krylov method 12073b09bd56SBarry Smith 12088cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 12098cf6037eSBarry Smith 121023067569SBarry Smith When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 121123067569SBarry Smith is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 121223067569SBarry Smith (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 121323067569SBarry Smith residual is computed at the end of each cycle. 121423067569SBarry Smith 12153b09bd56SBarry Smith Level: intermediate 12163b09bd56SBarry Smith 12178f87f92bSBarry Smith Concepts: multigrid/multilevel 12183b09bd56SBarry Smith 12198cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1220710315b6SLawrence Mitchell PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), 1221710315b6SLawrence Mitchell PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 122297177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 12230d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 12243b09bd56SBarry Smith M*/ 12253b09bd56SBarry Smith 12268cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 12274b9ad928SBarry Smith { 1228f3fbd535SBarry Smith PC_MG *mg; 1229f3fbd535SBarry Smith PetscErrorCode ierr; 1230f3fbd535SBarry Smith 12314b9ad928SBarry Smith PetscFunctionBegin; 1232b00a9115SJed Brown ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1233f3fbd535SBarry Smith pc->data = (void*)mg; 1234f3fbd535SBarry Smith mg->nlevels = -1; 123510eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 12362134b1e4SBarry Smith mg->galerkin = PC_MG_GALERKIN_NONE; 1237f3fbd535SBarry Smith 123837a44384SMark Adams pc->useAmat = PETSC_TRUE; 123937a44384SMark Adams 12404b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 12414b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1242a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 12434b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 12444b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 12454b9ad928SBarry Smith pc->ops->view = PCView_MG; 1246fb15c04eSBarry Smith 1247fb15c04eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 12484b9ad928SBarry Smith PetscFunctionReturn(0); 12494b9ad928SBarry Smith } 1250