1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 34b9ad928SBarry Smith /* 44b9ad928SBarry Smith Defines the multigrid preconditioner interface. 54b9ad928SBarry Smith */ 67c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 74b9ad928SBarry Smith 84b9ad928SBarry Smith 94b9ad928SBarry Smith #undef __FUNCT__ 109dcbbd2bSBarry Smith #define __FUNCT__ "PCMGMCycle_Private" 114d0a8057SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG **mglevels,PCRichardsonConvergedReason *reason) 124b9ad928SBarry Smith { 139dcbbd2bSBarry Smith PC_MG *mg = *mglevels,*mgc; 146849ba73SBarry Smith PetscErrorCode ierr; 150d353602SBarry Smith PetscInt cycles = (PetscInt) mg->cycles; 164b9ad928SBarry Smith 174b9ad928SBarry Smith PetscFunctionBegin; 184b9ad928SBarry Smith 1932cf1786SBarry Smith if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 200d353602SBarry Smith ierr = KSPSolve(mg->smoothd,mg->b,mg->x);CHKERRQ(ierr); /* pre-smooth */ 2132cf1786SBarry Smith if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 224b9ad928SBarry Smith if (mg->level) { /* not the coarsest grid */ 2332cf1786SBarry Smith if (mg->eventresidual) {ierr = PetscLogEventBegin(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);} 244b9ad928SBarry Smith ierr = (*mg->residual)(mg->A,mg->b,mg->x,mg->r);CHKERRQ(ierr); 2532cf1786SBarry Smith if (mg->eventresidual) {ierr = PetscLogEventEnd(mg->eventresidual,0,0,0,0);CHKERRQ(ierr);} 264b9ad928SBarry Smith 274b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 284d0a8057SBarry Smith if (mg->level == mg->levels-1 && mg->ttol && reason) { 294b9ad928SBarry Smith PetscReal rnorm; 304b9ad928SBarry Smith ierr = VecNorm(mg->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 434b9ad928SBarry Smith mgc = *(mglevels - 1); 4432cf1786SBarry Smith if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 454b9ad928SBarry Smith ierr = MatRestrict(mg->restrct,mg->r,mgc->b);CHKERRQ(ierr); 4632cf1786SBarry Smith if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 47efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 484b9ad928SBarry Smith while (cycles--) { 494d0a8057SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels-1,reason);CHKERRQ(ierr); 504b9ad928SBarry Smith } 5132cf1786SBarry Smith if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 524b9ad928SBarry Smith ierr = MatInterpolateAdd(mg->interpolate,mgc->x,mg->x,mg->x);CHKERRQ(ierr); 5332cf1786SBarry Smith if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5432cf1786SBarry Smith if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 550d353602SBarry Smith ierr = KSPSolve(mg->smoothu,mg->b,mg->x);CHKERRQ(ierr); /* post smooth */ 5632cf1786SBarry Smith if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 574b9ad928SBarry Smith } 584b9ad928SBarry Smith PetscFunctionReturn(0); 594b9ad928SBarry Smith } 604b9ad928SBarry Smith 614b9ad928SBarry Smith /* 629dcbbd2bSBarry Smith PCMGCreate_Private - Creates a PC_MG structure for use with the 634b9ad928SBarry Smith multigrid code. Level 0 is the coarsest. (But the 644b9ad928SBarry Smith finest level is stored first in the array). 654b9ad928SBarry Smith 664b9ad928SBarry Smith */ 674b9ad928SBarry Smith #undef __FUNCT__ 689dcbbd2bSBarry Smith #define __FUNCT__ "PCMGCreate_Private" 699dcbbd2bSBarry Smith static PetscErrorCode PCMGCreate_Private(MPI_Comm comm,PetscInt levels,PC pc,MPI_Comm *comms,PC_MG ***result) 704b9ad928SBarry Smith { 719dcbbd2bSBarry Smith PC_MG **mg; 726849ba73SBarry Smith PetscErrorCode ierr; 7379416396SBarry Smith PetscInt i; 7479416396SBarry Smith PetscMPIInt size; 75f69a0ea3SMatthew Knepley const char *prefix; 764b9ad928SBarry Smith PC ipc; 774b9ad928SBarry Smith 784b9ad928SBarry Smith PetscFunctionBegin; 799dcbbd2bSBarry Smith ierr = PetscMalloc(levels*sizeof(PC_MG*),&mg);CHKERRQ(ierr); 8038f2d2fdSLisandro Dalcin ierr = PetscLogObjectMemory(pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 814b9ad928SBarry Smith 824b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 834b9ad928SBarry Smith 844b9ad928SBarry Smith for (i=0; i<levels; i++) { 8538f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_MG,&mg[i]);CHKERRQ(ierr); 864b9ad928SBarry Smith mg[i]->level = i; 874b9ad928SBarry Smith mg[i]->levels = levels; 880d353602SBarry Smith mg[i]->cycles = PC_MG_CYCLE_V; 89c2be2410SBarry Smith mg[i]->galerkin = PETSC_FALSE; 90c2be2410SBarry Smith mg[i]->galerkinused = PETSC_FALSE; 9179416396SBarry Smith mg[i]->default_smoothu = 1; 9279416396SBarry Smith mg[i]->default_smoothd = 1; 934b9ad928SBarry Smith 944b9ad928SBarry Smith if (comms) comm = comms[i]; 954b9ad928SBarry Smith ierr = KSPCreate(comm,&mg[i]->smoothd);CHKERRQ(ierr); 961cee3971SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mg[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 9779416396SBarry Smith ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg[i]->default_smoothd);CHKERRQ(ierr); 984b9ad928SBarry Smith ierr = KSPSetOptionsPrefix(mg[i]->smoothd,prefix);CHKERRQ(ierr); 994b9ad928SBarry Smith 1004b9ad928SBarry Smith /* do special stuff for coarse grid */ 1014b9ad928SBarry Smith if (!i && levels > 1) { 1024b9ad928SBarry Smith ierr = KSPAppendOptionsPrefix(mg[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 1034b9ad928SBarry Smith 1044b9ad928SBarry Smith /* coarse solve is (redundant) LU by default */ 1054b9ad928SBarry Smith ierr = KSPSetType(mg[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 1064b9ad928SBarry Smith ierr = KSPGetPC(mg[0]->smoothd,&ipc);CHKERRQ(ierr); 1074b9ad928SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1084b9ad928SBarry Smith if (size > 1) { 1094b9ad928SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 1100d7810c8SBarry Smith } else { 1114b9ad928SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 1120d7810c8SBarry Smith } 1134b9ad928SBarry Smith 1144b9ad928SBarry Smith } else { 1152e70eadcSBarry Smith char tprefix[128]; 11613f74950SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 1172e70eadcSBarry Smith ierr = KSPAppendOptionsPrefix(mg[i]->smoothd,tprefix);CHKERRQ(ierr); 1184b9ad928SBarry Smith } 11952e6d16bSBarry Smith ierr = PetscLogObjectParent(pc,mg[i]->smoothd);CHKERRQ(ierr); 1204b9ad928SBarry Smith mg[i]->smoothu = mg[i]->smoothd; 1214b9ad928SBarry Smith mg[i]->rtol = 0.0; 12270441072SBarry Smith mg[i]->abstol = 0.0; 1234b9ad928SBarry Smith mg[i]->dtol = 0.0; 1244b9ad928SBarry Smith mg[i]->ttol = 0.0; 12532cf1786SBarry Smith mg[i]->eventsmoothsetup = 0; 12632cf1786SBarry Smith mg[i]->eventsmoothsolve = 0; 12732cf1786SBarry Smith mg[i]->eventresidual = 0; 12832cf1786SBarry Smith mg[i]->eventinterprestrict = 0; 1298cc2d5dfSBarry Smith mg[i]->cyclesperpcapply = 1; 1304b9ad928SBarry Smith } 1314b9ad928SBarry Smith *result = mg; 1324b9ad928SBarry Smith PetscFunctionReturn(0); 1334b9ad928SBarry Smith } 1344b9ad928SBarry Smith 1354b9ad928SBarry Smith #undef __FUNCT__ 136*14ca34e6SBarry Smith #define __FUNCT__ "PCMGSetSetup" 137*14ca34e6SBarry Smith /*@ 138*14ca34e6SBarry Smith PCMGSetSetup - sets a function that is called for each new setup of the MG PC that allows the user 139*14ca34e6SBarry Smith to provide a new interpolation/coarse matrix etc based on new solver matrix etc 140*14ca34e6SBarry Smith 141*14ca34e6SBarry Smith Collective on PC 142*14ca34e6SBarry Smith 143*14ca34e6SBarry Smith Input Parameters: 144*14ca34e6SBarry Smith + pc - the preconditioner context 145*14ca34e6SBarry Smith . setup - the function that is called each time PCSetUp() is called 146*14ca34e6SBarry Smith . setupdestroy - optional function that destroys the setup context 147*14ca34e6SBarry Smith - setupctx - optional context passed into the setup function 148*14ca34e6SBarry Smith 149*14ca34e6SBarry Smith Level: intermediate 150*14ca34e6SBarry Smith 151*14ca34e6SBarry Smith .keywords: PC, MG 152*14ca34e6SBarry Smith @*/ 153*14ca34e6SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetSetup(PC pc,PetscErrorCode (*setup)(PC,void*),PetscErrorCode (*setupdestroy)(PC,void*),void *setupctx) 154*14ca34e6SBarry Smith { 155*14ca34e6SBarry Smith PetscErrorCode ierr,(*f)(PC,PetscErrorCode (*setup)(PC,void*),PetscErrorCode (*setupdestroy)(PC,void*),void *setupctx); 156*14ca34e6SBarry Smith 157*14ca34e6SBarry Smith PetscFunctionBegin; 158*14ca34e6SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 159*14ca34e6SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCMGSetSetup_C",(void (**)(void))&f);CHKERRQ(ierr); 160*14ca34e6SBarry Smith if (f) { 161*14ca34e6SBarry Smith ierr = (*f)(pc,setup,setupdestroy,setupctx);CHKERRQ(ierr); 162*14ca34e6SBarry Smith } 163*14ca34e6SBarry Smith PetscFunctionReturn(0); 164*14ca34e6SBarry Smith } 165*14ca34e6SBarry Smith 166*14ca34e6SBarry Smith #undef __FUNCT__ 167*14ca34e6SBarry Smith #define __FUNCT__ "PCMGSetSetup_MG" 168*14ca34e6SBarry Smith PetscErrorCode PCMGSetSetup_MG(PC pc,PetscErrorCode (*setup)(PC,void*),PetscErrorCode (*setupdestroy)(PC,void*),void *setupctx) 169*14ca34e6SBarry Smith { 170*14ca34e6SBarry Smith PC_MG **mg = (PC_MG**)pc->data; 171*14ca34e6SBarry Smith 172*14ca34e6SBarry Smith PetscFunctionBegin; 173*14ca34e6SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 174*14ca34e6SBarry Smith mg[0]->setup = setup; 175*14ca34e6SBarry Smith mg[0]->setupdestroy = setupdestroy; 176*14ca34e6SBarry Smith mg[0]->setupctx = setupctx; 177*14ca34e6SBarry Smith PetscFunctionReturn(0); 178*14ca34e6SBarry Smith } 179*14ca34e6SBarry Smith 180*14ca34e6SBarry Smith #undef __FUNCT__ 1814b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_MG" 1826849ba73SBarry Smith static PetscErrorCode PCDestroy_MG(PC pc) 1834b9ad928SBarry Smith { 1849dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 1856849ba73SBarry Smith PetscErrorCode ierr; 18638f2d2fdSLisandro Dalcin PetscInt i,n; 1874b9ad928SBarry Smith 1884b9ad928SBarry Smith PetscFunctionBegin; 18938f2d2fdSLisandro Dalcin if (!mg) PetscFunctionReturn(0); 190*14ca34e6SBarry Smith 191*14ca34e6SBarry Smith if (mg[0]->setupdestroy) { 192*14ca34e6SBarry Smith ierr = (*mg[0]->setupdestroy)(pc,mg[0]->setupctx);CHKERRQ(ierr); 193*14ca34e6SBarry Smith } 194*14ca34e6SBarry Smith 19538f2d2fdSLisandro Dalcin n = mg[0]->levels; 196fccaa45eSBarry Smith for (i=0; i<n-1; i++) { 197630ba2eeSBarry Smith if (mg[i+1]->r) {ierr = VecDestroy(mg[i+1]->r);CHKERRQ(ierr);} 198fccaa45eSBarry Smith if (mg[i]->b) {ierr = VecDestroy(mg[i]->b);CHKERRQ(ierr);} 199fccaa45eSBarry Smith if (mg[i]->x) {ierr = VecDestroy(mg[i]->x);CHKERRQ(ierr);} 200fccaa45eSBarry Smith if (mg[i+1]->restrct) {ierr = MatDestroy(mg[i+1]->restrct);CHKERRQ(ierr);} 201fccaa45eSBarry Smith if (mg[i+1]->interpolate) {ierr = MatDestroy(mg[i+1]->interpolate);CHKERRQ(ierr);} 202fccaa45eSBarry Smith } 203fccaa45eSBarry Smith 2044b9ad928SBarry Smith for (i=0; i<n; i++) { 2054b9ad928SBarry Smith if (mg[i]->smoothd != mg[i]->smoothu) { 2064b9ad928SBarry Smith ierr = KSPDestroy(mg[i]->smoothd);CHKERRQ(ierr); 2074b9ad928SBarry Smith } 2084b9ad928SBarry Smith ierr = KSPDestroy(mg[i]->smoothu);CHKERRQ(ierr); 2094b9ad928SBarry Smith ierr = PetscFree(mg[i]);CHKERRQ(ierr); 2104b9ad928SBarry Smith } 2114b9ad928SBarry Smith ierr = PetscFree(mg);CHKERRQ(ierr); 2124b9ad928SBarry Smith PetscFunctionReturn(0); 2134b9ad928SBarry Smith } 2144b9ad928SBarry Smith 2154b9ad928SBarry Smith 2164b9ad928SBarry Smith 2179dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGACycle_Private(PC_MG**); 2181e2582c4SBarry Smith EXTERN PetscErrorCode PCMGFCycle_Private(PC,PC_MG**); 2199dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGKCycle_Private(PC_MG**); 2204b9ad928SBarry Smith 2214b9ad928SBarry Smith /* 2224b9ad928SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 2234b9ad928SBarry Smith or full cycle of multigrid. 2244b9ad928SBarry Smith 2254b9ad928SBarry Smith Note: 2269dcbbd2bSBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 2274b9ad928SBarry Smith */ 2284b9ad928SBarry Smith #undef __FUNCT__ 2294b9ad928SBarry Smith #define __FUNCT__ "PCApply_MG" 2306849ba73SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 2314b9ad928SBarry Smith { 2329dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 2336849ba73SBarry Smith PetscErrorCode ierr; 2348cc2d5dfSBarry Smith PetscInt levels = mg[0]->levels,i; 2354b9ad928SBarry Smith 2364b9ad928SBarry Smith PetscFunctionBegin; 2374b9ad928SBarry Smith mg[levels-1]->b = b; 2384b9ad928SBarry Smith mg[levels-1]->x = x; 2399dcbbd2bSBarry Smith if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 240efb30889SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 2418cc2d5dfSBarry Smith for (i=0; i<mg[0]->cyclesperpcapply; i++) { 2421e2582c4SBarry Smith ierr = PCMGMCycle_Private(pc,mg+levels-1,PETSC_NULL);CHKERRQ(ierr); 2434b9ad928SBarry Smith } 2448cc2d5dfSBarry Smith } 2459dcbbd2bSBarry Smith else if (mg[0]->am == PC_MG_ADDITIVE) { 2469dcbbd2bSBarry Smith ierr = PCMGACycle_Private(mg);CHKERRQ(ierr); 2474b9ad928SBarry Smith } 2489dcbbd2bSBarry Smith else if (mg[0]->am == PC_MG_KASKADE) { 2499dcbbd2bSBarry Smith ierr = PCMGKCycle_Private(mg);CHKERRQ(ierr); 2504b9ad928SBarry Smith } 2514b9ad928SBarry Smith else { 2521e2582c4SBarry Smith ierr = PCMGFCycle_Private(pc,mg);CHKERRQ(ierr); 2534b9ad928SBarry Smith } 2544b9ad928SBarry Smith PetscFunctionReturn(0); 2554b9ad928SBarry Smith } 2564b9ad928SBarry Smith 2574b9ad928SBarry Smith #undef __FUNCT__ 2584b9ad928SBarry Smith #define __FUNCT__ "PCApplyRichardson_MG" 2594d0a8057SBarry Smith static PetscErrorCode PCApplyRichardson_MG(PC pc,Vec b,Vec x,Vec w,PetscReal rtol,PetscReal abstol, PetscReal dtol,PetscInt its,PetscInt *outits,PCRichardsonConvergedReason *reason) 2604b9ad928SBarry Smith { 2619dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 262dfbe8321SBarry Smith PetscErrorCode ierr; 2634d0a8057SBarry Smith PetscInt levels = mg[0]->levels,i; 2644b9ad928SBarry Smith 2654b9ad928SBarry Smith PetscFunctionBegin; 2664b9ad928SBarry Smith mg[levels-1]->b = b; 2674b9ad928SBarry Smith mg[levels-1]->x = x; 2684b9ad928SBarry Smith 2694b9ad928SBarry Smith mg[levels-1]->rtol = rtol; 27070441072SBarry Smith mg[levels-1]->abstol = abstol; 2714b9ad928SBarry Smith mg[levels-1]->dtol = dtol; 2724b9ad928SBarry Smith if (rtol) { 2734b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 2744b9ad928SBarry Smith PetscReal rnorm; 2754b9ad928SBarry Smith ierr = (*mg[levels-1]->residual)(mg[levels-1]->A,b,x,w);CHKERRQ(ierr); 2764b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 27770441072SBarry Smith mg[levels-1]->ttol = PetscMax(rtol*rnorm,abstol); 27870441072SBarry Smith } else if (abstol) { 27970441072SBarry Smith mg[levels-1]->ttol = abstol; 2804b9ad928SBarry Smith } else { 2814b9ad928SBarry Smith mg[levels-1]->ttol = 0.0; 2824b9ad928SBarry Smith } 2834b9ad928SBarry Smith 2844d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 2854d0a8057SBarry Smith stop prematurely do to small residual */ 2864d0a8057SBarry Smith for (i=1; i<levels; i++) { 2874d0a8057SBarry Smith ierr = KSPSetTolerances(mg[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 2884d0a8057SBarry Smith if (mg[i]->smoothu != mg[i]->smoothd) { 2894d0a8057SBarry Smith ierr = KSPSetTolerances(mg[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 2904b9ad928SBarry Smith } 2914d0a8057SBarry Smith } 2924d0a8057SBarry Smith 2934d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 2944d0a8057SBarry Smith for (i=0; i<its; i++) { 2954d0a8057SBarry Smith ierr = PCMGMCycle_Private(pc,mg+levels-1,reason);CHKERRQ(ierr); 2964d0a8057SBarry Smith if (*reason) break; 2974d0a8057SBarry Smith } 2984d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 2994d0a8057SBarry Smith *outits = i; 3004b9ad928SBarry Smith PetscFunctionReturn(0); 3014b9ad928SBarry Smith } 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith #undef __FUNCT__ 3044b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_MG" 3056ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_MG(PC pc) 3064b9ad928SBarry Smith { 307dfbe8321SBarry Smith PetscErrorCode ierr; 3088cc2d5dfSBarry Smith PetscInt m,levels = 1,cycles; 3094b9ad928SBarry Smith PetscTruth flg; 3109dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 31190da95b0SMatthew Knepley PCMGType mgtype = PC_MG_ADDITIVE; 3121aa34eeaSBarry Smith PCMGCycleType mgctype; 3134b9ad928SBarry Smith 3144b9ad928SBarry Smith PetscFunctionBegin; 3154b9ad928SBarry Smith ierr = PetscOptionsHead("Multigrid options");CHKERRQ(ierr); 3164b9ad928SBarry Smith if (!pc->data) { 3179dcbbd2bSBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 31897177400SBarry Smith ierr = PCMGSetLevels(pc,levels,PETSC_NULL);CHKERRQ(ierr); 319cf502942SBarry Smith mg = (PC_MG**)pc->data; 3204b9ad928SBarry Smith } 321f2070a76SMatthew Knepley mgctype = (PCMGCycleType) mg[0]->cycles; 3220d353602SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 3234b9ad928SBarry Smith if (flg) { 3240d353602SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 3250d353602SBarry Smith }; 3269dcbbd2bSBarry Smith ierr = PetscOptionsName("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",&flg);CHKERRQ(ierr); 327c2be2410SBarry Smith if (flg) { 32897177400SBarry Smith ierr = PCMGSetGalerkin(pc);CHKERRQ(ierr); 329c2be2410SBarry Smith } 3309dcbbd2bSBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","PCMGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 3314b9ad928SBarry Smith if (flg) { 33297177400SBarry Smith ierr = PCMGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 3334b9ad928SBarry Smith } 3349dcbbd2bSBarry Smith ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","PCMGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 3354b9ad928SBarry Smith if (flg) { 33697177400SBarry Smith ierr = PCMGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 3374b9ad928SBarry Smith } 3389dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 3398cc2d5dfSBarry Smith if (flg) { 3408cc2d5dfSBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 3418cc2d5dfSBarry Smith } 3428cc2d5dfSBarry Smith if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 3438cc2d5dfSBarry Smith ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGSetLevels",mg[0]->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 3448cc2d5dfSBarry Smith if (flg) { 3458cc2d5dfSBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 3468cc2d5dfSBarry Smith } 3478cc2d5dfSBarry Smith } 3484b9ad928SBarry Smith ierr = PetscOptionsName("-pc_mg_log","Log times for each multigrid level","None",&flg);CHKERRQ(ierr); 3494b9ad928SBarry Smith if (flg) { 3504f5ab15aSBarry Smith PetscInt i; 3514b9ad928SBarry Smith char eventname[128]; 3524b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 3534b9ad928SBarry Smith levels = mg[0]->levels; 3544b9ad928SBarry Smith for (i=0; i<levels; i++) { 35532cf1786SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 3568cbcd9ccSBarry Smith ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsetup);CHKERRQ(ierr); 35732cf1786SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 3588cbcd9ccSBarry Smith ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventsmoothsolve);CHKERRQ(ierr); 35932cf1786SBarry Smith if (i) { 36032cf1786SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 361a3bc4eb9SBarry Smith ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventresidual);CHKERRQ(ierr); 36232cf1786SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 3638cbcd9ccSBarry Smith ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->cookie,&mg[i]->eventinterprestrict);CHKERRQ(ierr); 36432cf1786SBarry Smith } 3654b9ad928SBarry Smith } 3664b9ad928SBarry Smith } 3674b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 3684b9ad928SBarry Smith PetscFunctionReturn(0); 3694b9ad928SBarry Smith } 3704b9ad928SBarry Smith 3719dcbbd2bSBarry Smith const char *PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",0}; 3720d353602SBarry Smith const char *PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",0}; 3739dcbbd2bSBarry Smith 3744b9ad928SBarry Smith #undef __FUNCT__ 3754b9ad928SBarry Smith #define __FUNCT__ "PCView_MG" 3766849ba73SBarry Smith static PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 3774b9ad928SBarry Smith { 3789dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 379dfbe8321SBarry Smith PetscErrorCode ierr; 38079416396SBarry Smith PetscInt levels = mg[0]->levels,i; 38132077d6dSBarry Smith PetscTruth iascii; 3824b9ad928SBarry Smith 3834b9ad928SBarry Smith PetscFunctionBegin; 38432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 38532077d6dSBarry Smith if (iascii) { 386a8e6722dSMatthew Knepley ierr = PetscViewerASCIIPrintf(viewer," MG: type is %s, levels=%D cycles=%s\n", PCMGTypes[mg[0]->am],levels,(mg[0]->cycles == PC_MG_CYCLE_V) ? "v" : "w");CHKERRQ(ierr); 387d792f543SMatthew Knepley if (mg[0]->am == PC_MG_MULTIPLICATIVE) { 388d792f543SMatthew Knepley ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg[0]->cyclesperpcapply);CHKERRQ(ierr); 389d792f543SMatthew Knepley } 390c2be2410SBarry Smith if (mg[0]->galerkin) { 391c2be2410SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 392c2be2410SBarry Smith } 3934b9ad928SBarry Smith for (i=0; i<levels; i++) { 394b03c7568SBarry Smith if (!i) { 395a8e6722dSMatthew Knepley ierr = PetscViewerASCIIPrintf(viewer,"Coarse gride solver -- level %D presmooths=%D postsmooths=%D -----\n",i,mg[0]->default_smoothd,mg[0]->default_smoothu);CHKERRQ(ierr); 396b03c7568SBarry Smith } else { 397a8e6722dSMatthew Knepley ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothd);CHKERRQ(ierr); 398b03c7568SBarry Smith } 3994b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4004b9ad928SBarry Smith ierr = KSPView(mg[i]->smoothd,viewer);CHKERRQ(ierr); 4014b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 402b03c7568SBarry Smith if (i && mg[i]->smoothd == mg[i]->smoothu) { 4034b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 404b03c7568SBarry Smith } else if (i){ 405a8e6722dSMatthew Knepley ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D smooths=%D --------------------\n",i,mg[i]->default_smoothu);CHKERRQ(ierr); 4064b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4074b9ad928SBarry Smith ierr = KSPView(mg[i]->smoothu,viewer);CHKERRQ(ierr); 4084b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 4094b9ad928SBarry Smith } 4104b9ad928SBarry Smith } 4114b9ad928SBarry Smith } else { 41279a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PCMG",((PetscObject)viewer)->type_name); 4134b9ad928SBarry Smith } 4144b9ad928SBarry Smith PetscFunctionReturn(0); 4154b9ad928SBarry Smith } 4164b9ad928SBarry Smith 4174b9ad928SBarry Smith /* 4184b9ad928SBarry Smith Calls setup for the KSP on each level 4194b9ad928SBarry Smith */ 4204b9ad928SBarry Smith #undef __FUNCT__ 4214b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_MG" 4226849ba73SBarry Smith static PetscErrorCode PCSetUp_MG(PC pc) 4234b9ad928SBarry Smith { 4249dcbbd2bSBarry Smith PC_MG **mg = (PC_MG**)pc->data; 425dfbe8321SBarry Smith PetscErrorCode ierr; 42679416396SBarry Smith PetscInt i,n = mg[0]->levels; 42777122347SBarry Smith PC cpc,mpc; 428906ed7ccSBarry Smith PetscTruth preonly,lu,redundant,cholesky,monitor = PETSC_FALSE,dump,opsset; 42923d894e5SBarry Smith PetscViewerASCIIMonitor ascii; 43023d894e5SBarry Smith PetscViewer viewer = PETSC_NULL; 4314b9ad928SBarry Smith MPI_Comm comm; 432c2be2410SBarry Smith Mat dA,dB; 433c2be2410SBarry Smith MatStructure uflag; 4340a6bb862SBarry Smith Vec tvec; 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith PetscFunctionBegin; 437852665d3SBarry Smith 438*14ca34e6SBarry Smith if (mg[0]->setup) { 439*14ca34e6SBarry Smith ierr = (*mg[0]->setup)(pc,mg[0]->setupctx);CHKERRQ(ierr); 440*14ca34e6SBarry Smith } 441*14ca34e6SBarry Smith 44243fb2f97SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 44343fb2f97SBarry Smith /* so use those from global PC */ 44443fb2f97SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 445906ed7ccSBarry Smith ierr = KSPGetOperatorsSet(mg[n-1]->smoothd,PETSC_NULL,&opsset);CHKERRQ(ierr); 44643fb2f97SBarry Smith ierr = KSPGetPC(mg[0]->smoothd,&cpc);CHKERRQ(ierr); 44777122347SBarry Smith ierr = KSPGetPC(mg[n-1]->smoothd,&mpc);CHKERRQ(ierr); 44877122347SBarry Smith if (!opsset || ((cpc->setupcalled == 1) && (mpc->setupcalled == 2))) { 449852665d3SBarry 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); 4501cfe3bddSBarry Smith ierr = KSPSetOperators(mg[n-1]->smoothd,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 451852665d3SBarry Smith } 452852665d3SBarry Smith 453852665d3SBarry Smith if (mg[0]->galerkin) { 454852665d3SBarry Smith Mat B; 455852665d3SBarry Smith mg[0]->galerkinused = PETSC_TRUE; 456852665d3SBarry Smith /* currently only handle case where mat and pmat are the same on coarser levels */ 457852665d3SBarry Smith ierr = KSPGetOperators(mg[n-1]->smoothd,&dA,&dB,&uflag);CHKERRQ(ierr); 458852665d3SBarry Smith if (!pc->setupcalled) { 459852665d3SBarry Smith for (i=n-2; i>-1; i--) { 460852665d3SBarry Smith ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_INITIAL_MATRIX,1.0,&B);CHKERRQ(ierr); 461852665d3SBarry Smith ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 462906ed7ccSBarry Smith if (i != n-2) {ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr);} 463852665d3SBarry Smith dB = B; 464852665d3SBarry Smith } 465906ed7ccSBarry Smith ierr = PetscObjectDereference((PetscObject)dB);CHKERRQ(ierr); 466852665d3SBarry Smith } else { 467852665d3SBarry Smith for (i=n-2; i>-1; i--) { 468906ed7ccSBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&B,PETSC_NULL);CHKERRQ(ierr); 469852665d3SBarry Smith ierr = MatPtAP(dB,mg[i+1]->interpolate,MAT_REUSE_MATRIX,1.0,&B);CHKERRQ(ierr); 470852665d3SBarry Smith ierr = KSPSetOperators(mg[i]->smoothd,B,B,uflag);CHKERRQ(ierr); 471852665d3SBarry Smith dB = B; 472852665d3SBarry Smith } 473852665d3SBarry Smith } 474852665d3SBarry Smith } 475852665d3SBarry Smith 476958c9bccSBarry Smith if (!pc->setupcalled) { 4774b9ad928SBarry Smith ierr = PetscOptionsHasName(0,"-pc_mg_monitor",&monitor);CHKERRQ(ierr); 4784b9ad928SBarry Smith 479b03c7568SBarry Smith for (i=0; i<n; i++) { 4804b9ad928SBarry Smith if (monitor) { 4814b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothd,&comm);CHKERRQ(ierr); 48223d894e5SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 48323d894e5SBarry Smith ierr = KSPMonitorSet(mg[i]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 4844b9ad928SBarry Smith } 4854b9ad928SBarry Smith ierr = KSPSetFromOptions(mg[i]->smoothd);CHKERRQ(ierr); 4864b9ad928SBarry Smith } 487b03c7568SBarry Smith for (i=1; i<n; i++) { 488a98fc643SBarry Smith if (mg[i]->smoothu && (mg[i]->smoothu != mg[i]->smoothd)) { 4894b9ad928SBarry Smith if (monitor) { 4904b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[i]->smoothu,&comm);CHKERRQ(ierr); 49123d894e5SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n-i,&ascii);CHKERRQ(ierr); 49223d894e5SBarry Smith ierr = KSPMonitorSet(mg[i]->smoothu,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 4934b9ad928SBarry Smith } 4944b9ad928SBarry Smith ierr = KSPSetFromOptions(mg[i]->smoothu);CHKERRQ(ierr); 4954b9ad928SBarry Smith } 4964b9ad928SBarry Smith } 497fccaa45eSBarry Smith for (i=1; i<n; i++) { 4980cace4b0SBarry Smith if (!mg[i]->residual) { 4990cace4b0SBarry Smith Mat mat; 5000cace4b0SBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,PETSC_NULL,&mat,PETSC_NULL);CHKERRQ(ierr); 5010cace4b0SBarry Smith ierr = PCMGSetResidual(pc,i,PCMGDefaultResidual,mat);CHKERRQ(ierr); 5020cace4b0SBarry Smith } 503fccaa45eSBarry Smith if (mg[i]->restrct && !mg[i]->interpolate) { 504aea2a34eSBarry Smith ierr = PCMGSetInterpolation(pc,i,mg[i]->restrct);CHKERRQ(ierr); 505fccaa45eSBarry Smith } 506fccaa45eSBarry Smith if (!mg[i]->restrct && mg[i]->interpolate) { 50797177400SBarry Smith ierr = PCMGSetRestriction(pc,i,mg[i]->interpolate);CHKERRQ(ierr); 508fccaa45eSBarry Smith } 509fccaa45eSBarry Smith #if defined(PETSC_USE_DEBUG) 510fccaa45eSBarry Smith if (!mg[i]->restrct || !mg[i]->interpolate) { 511fccaa45eSBarry Smith SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Need to set restriction or interpolation on level %d",(int)i); 512fccaa45eSBarry Smith } 513fccaa45eSBarry Smith #endif 514fccaa45eSBarry Smith } 5150a6bb862SBarry Smith for (i=0; i<n-1; i++) { 51637b0e6c0SBarry Smith if (!mg[i]->b) { 517906ed7ccSBarry Smith Vec *vec; 518906ed7ccSBarry Smith ierr = KSPGetVecs(mg[i]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 519906ed7ccSBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 5204fdc791cSSatish Balay ierr = VecDestroy(*vec);CHKERRQ(ierr); 521906ed7ccSBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 52237b0e6c0SBarry Smith } 5236ca4d86aSHong Zhang if (!mg[i]->r && i) { 5240a6bb862SBarry Smith ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 52597177400SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 5260a6bb862SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 5270a6bb862SBarry Smith } 5280a6bb862SBarry Smith if (!mg[i]->x) { 5290a6bb862SBarry Smith ierr = VecDuplicate(mg[i]->b,&tvec);CHKERRQ(ierr); 53097177400SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 5310a6bb862SBarry Smith ierr = VecDestroy(tvec);CHKERRQ(ierr); 5320a6bb862SBarry Smith } 5330a6bb862SBarry Smith } 534dfef70bdSHong Zhang if (n != 1 && !mg[n-1]->r) { 535dfef70bdSHong Zhang /* PCMGSetR() on the finest level if user did not supply it */ 5360f77fa87SBarry Smith Vec *vec; 5370f77fa87SBarry Smith ierr = KSPGetVecs(mg[n-1]->smoothd,1,&vec,0,PETSC_NULL);CHKERRQ(ierr); 5380f77fa87SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 539797d2ccbSBarry Smith ierr = VecDestroy(*vec);CHKERRQ(ierr); 5400f77fa87SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 5410f77fa87SBarry Smith } 5424b9ad928SBarry Smith } 5434b9ad928SBarry Smith 544c2be2410SBarry Smith 5454b9ad928SBarry Smith for (i=1; i<n; i++) { 546b03c7568SBarry Smith if (mg[i]->smoothu == mg[i]->smoothd) { 547b03c7568SBarry Smith /* if doing only down then initial guess is zero */ 5484b9ad928SBarry Smith ierr = KSPSetInitialGuessNonzero(mg[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 549b03c7568SBarry Smith } 55032cf1786SBarry Smith if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 5514b9ad928SBarry Smith ierr = KSPSetUp(mg[i]->smoothd);CHKERRQ(ierr); 55232cf1786SBarry Smith if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 5534b9ad928SBarry Smith } 554b03c7568SBarry Smith for (i=1; i<n; i++) { 5554b9ad928SBarry Smith if (mg[i]->smoothu && mg[i]->smoothu != mg[i]->smoothd) { 556906ed7ccSBarry Smith Mat downmat,downpmat; 55797f1f81fSBarry Smith MatStructure matflag; 558906ed7ccSBarry Smith PetscTruth opsset; 55997f1f81fSBarry Smith 56097f1f81fSBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 561906ed7ccSBarry Smith ierr = KSPGetOperatorsSet(mg[i]->smoothu,&opsset,PETSC_NULL);CHKERRQ(ierr); 562906ed7ccSBarry Smith if (!opsset) { 563906ed7ccSBarry Smith ierr = KSPGetOperators(mg[i]->smoothd,&downmat,&downpmat,&matflag);CHKERRQ(ierr); 56497f1f81fSBarry Smith ierr = KSPSetOperators(mg[i]->smoothu,downmat,downpmat,matflag);CHKERRQ(ierr); 56597f1f81fSBarry Smith } 56697f1f81fSBarry Smith 5674b9ad928SBarry Smith ierr = KSPSetInitialGuessNonzero(mg[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 56832cf1786SBarry Smith if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 5694b9ad928SBarry Smith ierr = KSPSetUp(mg[i]->smoothu);CHKERRQ(ierr); 57032cf1786SBarry Smith if (mg[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 5714b9ad928SBarry Smith } 5724b9ad928SBarry Smith } 5734b9ad928SBarry Smith 5744b9ad928SBarry Smith /* 5754b9ad928SBarry Smith If coarse solver is not direct method then DO NOT USE preonly 5764b9ad928SBarry Smith */ 5774b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)mg[0]->smoothd,KSPPREONLY,&preonly);CHKERRQ(ierr); 5784b9ad928SBarry Smith if (preonly) { 5794b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCLU,&lu);CHKERRQ(ierr); 5804b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCREDUNDANT,&redundant);CHKERRQ(ierr); 58168eff7e6SBarry Smith ierr = PetscTypeCompare((PetscObject)cpc,PCCHOLESKY,&cholesky);CHKERRQ(ierr); 58268eff7e6SBarry Smith if (!lu && !redundant && !cholesky) { 5834b9ad928SBarry Smith ierr = KSPSetType(mg[0]->smoothd,KSPGMRES);CHKERRQ(ierr); 5844b9ad928SBarry Smith } 5854b9ad928SBarry Smith } 5864b9ad928SBarry Smith 587958c9bccSBarry Smith if (!pc->setupcalled) { 5884b9ad928SBarry Smith if (monitor) { 5894b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)mg[0]->smoothd,&comm);CHKERRQ(ierr); 59023d894e5SBarry Smith ierr = PetscViewerASCIIMonitorCreate(comm,"stdout",n,&ascii);CHKERRQ(ierr); 59123d894e5SBarry Smith ierr = KSPMonitorSet(mg[0]->smoothd,KSPMonitorDefault,ascii,(PetscErrorCode(*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 5924b9ad928SBarry Smith } 5934b9ad928SBarry Smith ierr = KSPSetFromOptions(mg[0]->smoothd);CHKERRQ(ierr); 5944b9ad928SBarry Smith } 5954b9ad928SBarry Smith 59632cf1786SBarry Smith if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 5974b9ad928SBarry Smith ierr = KSPSetUp(mg[0]->smoothd);CHKERRQ(ierr); 59832cf1786SBarry Smith if (mg[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mg[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 5994b9ad928SBarry Smith 6004b9ad928SBarry Smith /* 6016805f65bSBarry Smith Dump the interpolation/restriction matrices plus the 6024b9ad928SBarry Smith Jacobian/stiffness on each level. This allows Matlab users to 6036805f65bSBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 6046805f65bSBarry Smith 6056805f65bSBarry Smith Only support one or the other at the same time. 6066805f65bSBarry Smith */ 6076805f65bSBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 6087adad957SLisandro Dalcin ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump);CHKERRQ(ierr); 6094b9ad928SBarry Smith if (dump) { 6107adad957SLisandro Dalcin viewer = PETSC_VIEWER_SOCKET_(((PetscObject)pc)->comm); 6114b9ad928SBarry Smith } 612c45a1595SBarry Smith #endif 6137adad957SLisandro Dalcin ierr = PetscOptionsHasName(((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump);CHKERRQ(ierr); 614c2be2410SBarry Smith if (dump) { 6157adad957SLisandro Dalcin viewer = PETSC_VIEWER_BINARY_(((PetscObject)pc)->comm); 6166805f65bSBarry Smith } 6176805f65bSBarry Smith 6186805f65bSBarry Smith if (viewer) { 619c2be2410SBarry Smith for (i=1; i<n; i++) { 6206805f65bSBarry Smith ierr = MatView(mg[i]->restrct,viewer);CHKERRQ(ierr); 621c2be2410SBarry Smith } 622c2be2410SBarry Smith for (i=0; i<n; i++) { 623c2be2410SBarry Smith ierr = KSPGetPC(mg[i]->smoothd,&pc);CHKERRQ(ierr); 6246805f65bSBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 625c2be2410SBarry Smith } 626c2be2410SBarry Smith } 6274b9ad928SBarry Smith PetscFunctionReturn(0); 6284b9ad928SBarry Smith } 6294b9ad928SBarry Smith 6304b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 6314b9ad928SBarry Smith 6324b9ad928SBarry Smith #undef __FUNCT__ 6339dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetLevels" 6344b9ad928SBarry Smith /*@C 63597177400SBarry Smith PCMGSetLevels - Sets the number of levels to use with MG. 6364b9ad928SBarry Smith Must be called before any other MG routine. 6374b9ad928SBarry Smith 6384b9ad928SBarry Smith Collective on PC 6394b9ad928SBarry Smith 6404b9ad928SBarry Smith Input Parameters: 6414b9ad928SBarry Smith + pc - the preconditioner context 6424b9ad928SBarry Smith . levels - the number of levels 6434b9ad928SBarry Smith - comms - optional communicators for each level; this is to allow solving the coarser problems 6444b9ad928SBarry Smith on smaller sets of processors. Use PETSC_NULL_OBJECT for default in Fortran 6454b9ad928SBarry Smith 6464b9ad928SBarry Smith Level: intermediate 6474b9ad928SBarry Smith 6484b9ad928SBarry Smith Notes: 6494b9ad928SBarry Smith If the number of levels is one then the multigrid uses the -mg_levels prefix 6504b9ad928SBarry Smith for setting the level options rather than the -mg_coarse prefix. 6514b9ad928SBarry Smith 6524b9ad928SBarry Smith .keywords: MG, set, levels, multigrid 6534b9ad928SBarry Smith 65497177400SBarry Smith .seealso: PCMGSetType(), PCMGGetLevels() 6554b9ad928SBarry Smith @*/ 65697177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 6574b9ad928SBarry Smith { 658dfbe8321SBarry Smith PetscErrorCode ierr; 659ada7143aSSatish Balay PC_MG **mg=0; 6604b9ad928SBarry Smith 6614b9ad928SBarry Smith PetscFunctionBegin; 6624482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 6634b9ad928SBarry Smith 6644b9ad928SBarry Smith if (pc->data) { 6651302d50aSBarry Smith SETERRQ(PETSC_ERR_ORDER,"Number levels already set for MG\n\ 66697177400SBarry Smith make sure that you call PCMGSetLevels() before KSPSetFromOptions()"); 6674b9ad928SBarry Smith } 6687adad957SLisandro Dalcin ierr = PCMGCreate_Private(((PetscObject)pc)->comm,levels,pc,comms,&mg);CHKERRQ(ierr); 6699dcbbd2bSBarry Smith mg[0]->am = PC_MG_MULTIPLICATIVE; 6704b9ad928SBarry Smith pc->data = (void*)mg; 6714b9ad928SBarry Smith pc->ops->applyrichardson = PCApplyRichardson_MG; 6724b9ad928SBarry Smith PetscFunctionReturn(0); 6734b9ad928SBarry Smith } 6744b9ad928SBarry Smith 6754b9ad928SBarry Smith #undef __FUNCT__ 6769dcbbd2bSBarry Smith #define __FUNCT__ "PCMGGetLevels" 6774b9ad928SBarry Smith /*@ 67897177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 6794b9ad928SBarry Smith 6804b9ad928SBarry Smith Not Collective 6814b9ad928SBarry Smith 6824b9ad928SBarry Smith Input Parameter: 6834b9ad928SBarry Smith . pc - the preconditioner context 6844b9ad928SBarry Smith 6854b9ad928SBarry Smith Output parameter: 6864b9ad928SBarry Smith . levels - the number of levels 6874b9ad928SBarry Smith 6884b9ad928SBarry Smith Level: advanced 6894b9ad928SBarry Smith 6904b9ad928SBarry Smith .keywords: MG, get, levels, multigrid 6914b9ad928SBarry Smith 69297177400SBarry Smith .seealso: PCMGSetLevels() 6934b9ad928SBarry Smith @*/ 69497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetLevels(PC pc,PetscInt *levels) 6954b9ad928SBarry Smith { 6969dcbbd2bSBarry Smith PC_MG **mg; 6974b9ad928SBarry Smith 6984b9ad928SBarry Smith PetscFunctionBegin; 6994482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 7004482741eSBarry Smith PetscValidIntPointer(levels,2); 7014b9ad928SBarry Smith 7029dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 7034b9ad928SBarry Smith *levels = mg[0]->levels; 7044b9ad928SBarry Smith PetscFunctionReturn(0); 7054b9ad928SBarry Smith } 7064b9ad928SBarry Smith 7074b9ad928SBarry Smith #undef __FUNCT__ 7089dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetType" 7094b9ad928SBarry Smith /*@ 71097177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 7114b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 7124b9ad928SBarry Smith 7134b9ad928SBarry Smith Collective on PC 7144b9ad928SBarry Smith 7154b9ad928SBarry Smith Input Parameters: 7164b9ad928SBarry Smith + pc - the preconditioner context 7179dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 7189dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 7194b9ad928SBarry Smith 7204b9ad928SBarry Smith Options Database Key: 7214b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 7224b9ad928SBarry Smith additive, full, kaskade 7234b9ad928SBarry Smith 7244b9ad928SBarry Smith Level: advanced 7254b9ad928SBarry Smith 7264b9ad928SBarry Smith .keywords: MG, set, method, multiplicative, additive, full, Kaskade, multigrid 7274b9ad928SBarry Smith 72897177400SBarry Smith .seealso: PCMGSetLevels() 7294b9ad928SBarry Smith @*/ 7309dcbbd2bSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetType(PC pc,PCMGType form) 7314b9ad928SBarry Smith { 7329dcbbd2bSBarry Smith PC_MG **mg; 7334b9ad928SBarry Smith 7344b9ad928SBarry Smith PetscFunctionBegin; 7354482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 7369dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 7374b9ad928SBarry Smith 7384b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 7394b9ad928SBarry Smith mg[0]->am = form; 7409dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 7414b9ad928SBarry Smith else pc->ops->applyrichardson = 0; 7424b9ad928SBarry Smith PetscFunctionReturn(0); 7434b9ad928SBarry Smith } 7444b9ad928SBarry Smith 7454b9ad928SBarry Smith #undef __FUNCT__ 7460d353602SBarry Smith #define __FUNCT__ "PCMGSetCycleType" 7474b9ad928SBarry Smith /*@ 7480d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 7494b9ad928SBarry Smith complicated cycling. 7504b9ad928SBarry Smith 7514b9ad928SBarry Smith Collective on PC 7524b9ad928SBarry Smith 7534b9ad928SBarry Smith Input Parameters: 754c2be2410SBarry Smith + pc - the multigrid context 7550d353602SBarry Smith - PC_MG_CYCLE_V or PC_MG_CYCLE_W 7564b9ad928SBarry Smith 7574b9ad928SBarry Smith Options Database Key: 7580d353602SBarry Smith $ -pc_mg_cycle_type v or w 7594b9ad928SBarry Smith 7604b9ad928SBarry Smith Level: advanced 7614b9ad928SBarry Smith 7624b9ad928SBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 7634b9ad928SBarry Smith 7640d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 7654b9ad928SBarry Smith @*/ 7660d353602SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetCycleType(PC pc,PCMGCycleType n) 7674b9ad928SBarry Smith { 7689dcbbd2bSBarry Smith PC_MG **mg; 76979416396SBarry Smith PetscInt i,levels; 7704b9ad928SBarry Smith 7714b9ad928SBarry Smith PetscFunctionBegin; 7724482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 7739dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 7744b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 7754b9ad928SBarry Smith levels = mg[0]->levels; 7764b9ad928SBarry Smith 7774b9ad928SBarry Smith for (i=0; i<levels; i++) { 7784b9ad928SBarry Smith mg[i]->cycles = n; 7794b9ad928SBarry Smith } 7804b9ad928SBarry Smith PetscFunctionReturn(0); 7814b9ad928SBarry Smith } 7824b9ad928SBarry Smith 7834b9ad928SBarry Smith #undef __FUNCT__ 7848cc2d5dfSBarry Smith #define __FUNCT__ "PCMGMultiplicativeSetCycles" 7858cc2d5dfSBarry Smith /*@ 7868cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 7878cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 7888cc2d5dfSBarry Smith 7898cc2d5dfSBarry Smith Collective on PC 7908cc2d5dfSBarry Smith 7918cc2d5dfSBarry Smith Input Parameters: 7928cc2d5dfSBarry Smith + pc - the multigrid context 7938cc2d5dfSBarry Smith - n - number of cycles (default is 1) 7948cc2d5dfSBarry Smith 7958cc2d5dfSBarry Smith Options Database Key: 7968cc2d5dfSBarry Smith $ -pc_mg_multiplicative_cycles n 7978cc2d5dfSBarry Smith 7988cc2d5dfSBarry Smith Level: advanced 7998cc2d5dfSBarry Smith 8008cc2d5dfSBarry Smith Notes: This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 8018cc2d5dfSBarry Smith 8028cc2d5dfSBarry Smith .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 8038cc2d5dfSBarry Smith 8048cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 8058cc2d5dfSBarry Smith @*/ 8068cc2d5dfSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 8078cc2d5dfSBarry Smith { 8088cc2d5dfSBarry Smith PC_MG **mg; 8098cc2d5dfSBarry Smith PetscInt i,levels; 8108cc2d5dfSBarry Smith 8118cc2d5dfSBarry Smith PetscFunctionBegin; 8128cc2d5dfSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8138cc2d5dfSBarry Smith mg = (PC_MG**)pc->data; 8148cc2d5dfSBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 8158cc2d5dfSBarry Smith levels = mg[0]->levels; 8168cc2d5dfSBarry Smith 8178cc2d5dfSBarry Smith for (i=0; i<levels; i++) { 8188cc2d5dfSBarry Smith mg[i]->cyclesperpcapply = n; 8198cc2d5dfSBarry Smith } 8208cc2d5dfSBarry Smith PetscFunctionReturn(0); 8218cc2d5dfSBarry Smith } 8228cc2d5dfSBarry Smith 8238cc2d5dfSBarry Smith #undef __FUNCT__ 8249dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetGalerkin" 825c2be2410SBarry Smith /*@ 82697177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 827c2be2410SBarry Smith finest grid via the Galerkin process: A_i-1 = r_i * A_i * r_i^t 828c2be2410SBarry Smith 829c2be2410SBarry Smith Collective on PC 830c2be2410SBarry Smith 831c2be2410SBarry Smith Input Parameters: 8323fc8bf9cSBarry Smith . pc - the multigrid context 833c2be2410SBarry Smith 834c2be2410SBarry Smith Options Database Key: 835c2be2410SBarry Smith $ -pc_mg_galerkin 836c2be2410SBarry Smith 837c2be2410SBarry Smith Level: intermediate 838c2be2410SBarry Smith 839c2be2410SBarry Smith .keywords: MG, set, Galerkin 840c2be2410SBarry Smith 8413fc8bf9cSBarry Smith .seealso: PCMGGetGalerkin() 8423fc8bf9cSBarry Smith 843c2be2410SBarry Smith @*/ 84497177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetGalerkin(PC pc) 845c2be2410SBarry Smith { 8469dcbbd2bSBarry Smith PC_MG **mg; 847c2be2410SBarry Smith PetscInt i,levels; 848c2be2410SBarry Smith 849c2be2410SBarry Smith PetscFunctionBegin; 850c2be2410SBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8519dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 852c2be2410SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 853c2be2410SBarry Smith levels = mg[0]->levels; 854c2be2410SBarry Smith 855c2be2410SBarry Smith for (i=0; i<levels; i++) { 856c2be2410SBarry Smith mg[i]->galerkin = PETSC_TRUE; 857c2be2410SBarry Smith } 858c2be2410SBarry Smith PetscFunctionReturn(0); 859c2be2410SBarry Smith } 860c2be2410SBarry Smith 861c2be2410SBarry Smith #undef __FUNCT__ 8623fc8bf9cSBarry Smith #define __FUNCT__ "PCMGGetGalerkin" 8633fc8bf9cSBarry Smith /*@ 8643fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 8653fc8bf9cSBarry Smith A_i-1 = r_i * A_i * r_i^t 8663fc8bf9cSBarry Smith 8673fc8bf9cSBarry Smith Not Collective 8683fc8bf9cSBarry Smith 8693fc8bf9cSBarry Smith Input Parameter: 8703fc8bf9cSBarry Smith . pc - the multigrid context 8713fc8bf9cSBarry Smith 8723fc8bf9cSBarry Smith Output Parameter: 8733fc8bf9cSBarry Smith . gelerkin - PETSC_TRUE or PETSC_FALSE 8743fc8bf9cSBarry Smith 8753fc8bf9cSBarry Smith Options Database Key: 8763fc8bf9cSBarry Smith $ -pc_mg_galerkin 8773fc8bf9cSBarry Smith 8783fc8bf9cSBarry Smith Level: intermediate 8793fc8bf9cSBarry Smith 8803fc8bf9cSBarry Smith .keywords: MG, set, Galerkin 8813fc8bf9cSBarry Smith 8823fc8bf9cSBarry Smith .seealso: PCMGSetGalerkin() 8833fc8bf9cSBarry Smith 8843fc8bf9cSBarry Smith @*/ 8853fc8bf9cSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGGetGalerkin(PC pc,PetscTruth *galerkin) 8863fc8bf9cSBarry Smith { 8873fc8bf9cSBarry Smith PC_MG **mg; 8883fc8bf9cSBarry Smith 8893fc8bf9cSBarry Smith PetscFunctionBegin; 8903fc8bf9cSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 8913fc8bf9cSBarry Smith mg = (PC_MG**)pc->data; 8923fc8bf9cSBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 8933fc8bf9cSBarry Smith *galerkin = mg[0]->galerkin; 8943fc8bf9cSBarry Smith PetscFunctionReturn(0); 8953fc8bf9cSBarry Smith } 8963fc8bf9cSBarry Smith 8973fc8bf9cSBarry Smith #undef __FUNCT__ 8989dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothDown" 8994b9ad928SBarry Smith /*@ 90097177400SBarry Smith PCMGSetNumberSmoothDown - Sets the number of pre-smoothing steps to 90197177400SBarry Smith use on all levels. Use PCMGGetSmootherDown() to set different 9024b9ad928SBarry Smith pre-smoothing steps on different levels. 9034b9ad928SBarry Smith 9044b9ad928SBarry Smith Collective on PC 9054b9ad928SBarry Smith 9064b9ad928SBarry Smith Input Parameters: 9074b9ad928SBarry Smith + mg - the multigrid context 9084b9ad928SBarry Smith - n - the number of smoothing steps 9094b9ad928SBarry Smith 9104b9ad928SBarry Smith Options Database Key: 9114b9ad928SBarry Smith . -pc_mg_smoothdown <n> - Sets number of pre-smoothing steps 9124b9ad928SBarry Smith 9134b9ad928SBarry Smith Level: advanced 9144b9ad928SBarry Smith 9154b9ad928SBarry Smith .keywords: MG, smooth, down, pre-smoothing, steps, multigrid 9164b9ad928SBarry Smith 91797177400SBarry Smith .seealso: PCMGSetNumberSmoothUp() 9184b9ad928SBarry Smith @*/ 91997177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothDown(PC pc,PetscInt n) 9204b9ad928SBarry Smith { 9219dcbbd2bSBarry Smith PC_MG **mg; 9226849ba73SBarry Smith PetscErrorCode ierr; 92379416396SBarry Smith PetscInt i,levels; 9244b9ad928SBarry Smith 9254b9ad928SBarry Smith PetscFunctionBegin; 9264482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 9279dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 9284b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 9294b9ad928SBarry Smith levels = mg[0]->levels; 9304b9ad928SBarry Smith 931b05257ddSBarry Smith for (i=1; i<levels; i++) { 9324b9ad928SBarry Smith /* make sure smoother up and down are different */ 93397177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 9344b9ad928SBarry Smith ierr = KSPSetTolerances(mg[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 9354b9ad928SBarry Smith mg[i]->default_smoothd = n; 9364b9ad928SBarry Smith } 9374b9ad928SBarry Smith PetscFunctionReturn(0); 9384b9ad928SBarry Smith } 9394b9ad928SBarry Smith 9404b9ad928SBarry Smith #undef __FUNCT__ 9419dcbbd2bSBarry Smith #define __FUNCT__ "PCMGSetNumberSmoothUp" 9424b9ad928SBarry Smith /*@ 94397177400SBarry Smith PCMGSetNumberSmoothUp - Sets the number of post-smoothing steps to use 94497177400SBarry Smith on all levels. Use PCMGGetSmootherUp() to set different numbers of 9454b9ad928SBarry Smith post-smoothing steps on different levels. 9464b9ad928SBarry Smith 9474b9ad928SBarry Smith Collective on PC 9484b9ad928SBarry Smith 9494b9ad928SBarry Smith Input Parameters: 9504b9ad928SBarry Smith + mg - the multigrid context 9514b9ad928SBarry Smith - n - the number of smoothing steps 9524b9ad928SBarry Smith 9534b9ad928SBarry Smith Options Database Key: 9544b9ad928SBarry Smith . -pc_mg_smoothup <n> - Sets number of post-smoothing steps 9554b9ad928SBarry Smith 9564b9ad928SBarry Smith Level: advanced 9574b9ad928SBarry Smith 9584b9ad928SBarry Smith Note: this does not set a value on the coarsest grid, since we assume that 959a8c7a070SBarry Smith there is no separate smooth up on the coarsest grid. 9604b9ad928SBarry Smith 9614b9ad928SBarry Smith .keywords: MG, smooth, up, post-smoothing, steps, multigrid 9624b9ad928SBarry Smith 96397177400SBarry Smith .seealso: PCMGSetNumberSmoothDown() 9644b9ad928SBarry Smith @*/ 96597177400SBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCMGSetNumberSmoothUp(PC pc,PetscInt n) 9664b9ad928SBarry Smith { 9679dcbbd2bSBarry Smith PC_MG **mg; 9686849ba73SBarry Smith PetscErrorCode ierr; 96979416396SBarry Smith PetscInt i,levels; 9704b9ad928SBarry Smith 9714b9ad928SBarry Smith PetscFunctionBegin; 9724482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 9739dcbbd2bSBarry Smith mg = (PC_MG**)pc->data; 9744b9ad928SBarry Smith if (!mg) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 9754b9ad928SBarry Smith levels = mg[0]->levels; 9764b9ad928SBarry Smith 9774b9ad928SBarry Smith for (i=1; i<levels; i++) { 9784b9ad928SBarry Smith /* make sure smoother up and down are different */ 97997177400SBarry Smith ierr = PCMGGetSmootherUp(pc,i,PETSC_NULL);CHKERRQ(ierr); 9804b9ad928SBarry Smith ierr = KSPSetTolerances(mg[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 9814b9ad928SBarry Smith mg[i]->default_smoothu = n; 9824b9ad928SBarry Smith } 9834b9ad928SBarry Smith PetscFunctionReturn(0); 9844b9ad928SBarry Smith } 9854b9ad928SBarry Smith 9864b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 9874b9ad928SBarry Smith 9883b09bd56SBarry Smith /*MC 989ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 9903b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 9913b09bd56SBarry Smith 9923b09bd56SBarry Smith Options Database Keys: 9933b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 9940d353602SBarry Smith . -pc_mg_cycles v or w 99579416396SBarry Smith . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 9963b09bd56SBarry Smith . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 9973b09bd56SBarry Smith . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 9983b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 9993b09bd56SBarry Smith . -pc_mg_monitor - print information on the multigrid convergence 100068eff7e6SBarry Smith . -pc_mg_galerkin - use Galerkin process to compute coarser operators 10013b09bd56SBarry Smith - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 10023b09bd56SBarry Smith to the Socket viewer for reading from Matlab. 10033b09bd56SBarry Smith 10043b09bd56SBarry Smith Notes: 10053b09bd56SBarry Smith 10063b09bd56SBarry Smith Level: intermediate 10073b09bd56SBarry Smith 10088f87f92bSBarry Smith Concepts: multigrid/multilevel 10093b09bd56SBarry Smith 10103b09bd56SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 10110d353602SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), PCMGSetNumberSmoothDown(), 101297177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 101397177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 10140d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 10153b09bd56SBarry Smith M*/ 10163b09bd56SBarry Smith 10174b9ad928SBarry Smith EXTERN_C_BEGIN 10184b9ad928SBarry Smith #undef __FUNCT__ 10194b9ad928SBarry Smith #define __FUNCT__ "PCCreate_MG" 1020dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_MG(PC pc) 10214b9ad928SBarry Smith { 1022*14ca34e6SBarry Smith PetscErrorCode ierr; 1023*14ca34e6SBarry Smith 10244b9ad928SBarry Smith PetscFunctionBegin; 10254b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 10264b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 10274b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 10284b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 10294b9ad928SBarry Smith pc->ops->view = PCView_MG; 10304b9ad928SBarry Smith 1031*14ca34e6SBarry Smith 1032*14ca34e6SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCMGSetSetup_C","PCMGSetSetup_MG",PCMGSetSetup_MG);CHKERRQ(ierr); 10334b9ad928SBarry Smith pc->data = (void*)0; 10344b9ad928SBarry Smith PetscFunctionReturn(0); 10354b9ad928SBarry Smith } 10364b9ad928SBarry Smith EXTERN_C_END 1037