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> 7391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*); 84b9ad928SBarry Smith 9*f3b08a26SMatthew G. Knepley /* 10*f3b08a26SMatthew G. Knepley Contains the list of registered coarse space construction routines 11*f3b08a26SMatthew G. Knepley */ 12*f3b08a26SMatthew G. Knepley PetscFunctionList PCMGCoarseList = NULL; 13*f3b08a26SMatthew G. Knepley 1431567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 154b9ad928SBarry Smith { 1631567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1731567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 186849ba73SBarry Smith PetscErrorCode ierr; 1931567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 204b9ad928SBarry Smith 214b9ad928SBarry Smith PetscFunctionBegin; 2263e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2331567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 24c0decd05SBarry Smith ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr); 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 */ 60c0decd05SBarry Smith ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr); 6163e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 624b9ad928SBarry Smith } 634b9ad928SBarry Smith PetscFunctionReturn(0); 644b9ad928SBarry Smith } 654b9ad928SBarry Smith 66ace3abfcSBarry 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) 674b9ad928SBarry Smith { 68f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 69f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 70dfbe8321SBarry Smith PetscErrorCode ierr; 71391689abSStefano Zampini PC tpc; 72391689abSStefano Zampini PetscBool changeu,changed; 73f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 744b9ad928SBarry Smith 754b9ad928SBarry Smith PetscFunctionBegin; 76a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 77a762d673SBarry Smith for (i=0; i<levels; i++) { 78a762d673SBarry Smith if (!mglevels[i]->A) { 79a762d673SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 80a762d673SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 81a762d673SBarry Smith } 82a762d673SBarry Smith } 83391689abSStefano Zampini 84391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 85391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 86391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 87391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 88391689abSStefano Zampini if (!changed && !changeu) { 89391689abSStefano Zampini ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 90f3fbd535SBarry Smith mglevels[levels-1]->b = b; 91391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 92391689abSStefano Zampini if (!mglevels[levels-1]->b) { 93391689abSStefano Zampini Vec *vec; 94391689abSStefano Zampini 95391689abSStefano Zampini ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 96391689abSStefano Zampini mglevels[levels-1]->b = *vec; 977ae21d43SStefano Zampini ierr = PetscFree(vec);CHKERRQ(ierr); 98391689abSStefano Zampini } 99391689abSStefano Zampini ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 100391689abSStefano Zampini } 101f3fbd535SBarry Smith mglevels[levels-1]->x = x; 1024b9ad928SBarry Smith 10331567311SBarry Smith mg->rtol = rtol; 10431567311SBarry Smith mg->abstol = abstol; 10531567311SBarry Smith mg->dtol = dtol; 1064b9ad928SBarry Smith if (rtol) { 1074b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 1084b9ad928SBarry Smith PetscReal rnorm; 1097319c654SBarry Smith if (zeroguess) { 1107319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 1117319c654SBarry Smith } else { 112f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 1134b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 1147319c654SBarry Smith } 11531567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 1162fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 1172fa5cd67SKarl Rupp else mg->ttol = 0.0; 1184b9ad928SBarry Smith 1194d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 120336babb1SJed Brown stop prematurely due to small residual */ 1214d0a8057SBarry Smith for (i=1; i<levels; i++) { 122f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 123f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 12423067569SBarry Smith /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 12523067569SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 126f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 1274b9ad928SBarry Smith } 1284d0a8057SBarry Smith } 1294d0a8057SBarry Smith 1304d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1314d0a8057SBarry Smith for (i=0; i<its; i++) { 132f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1334d0a8057SBarry Smith if (*reason) break; 1344d0a8057SBarry Smith } 1354d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1364d0a8057SBarry Smith *outits = i; 137391689abSStefano Zampini if (!changed && !changeu) mglevels[levels-1]->b = NULL; 1384b9ad928SBarry Smith PetscFunctionReturn(0); 1394b9ad928SBarry Smith } 1404b9ad928SBarry Smith 1413aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 1423aeaf226SBarry Smith { 1433aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1443aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1453aeaf226SBarry Smith PetscErrorCode ierr; 146*f3b08a26SMatthew G. Knepley PetscInt i,c,n; 1473aeaf226SBarry Smith 1483aeaf226SBarry Smith PetscFunctionBegin; 1493aeaf226SBarry Smith if (mglevels) { 1503aeaf226SBarry Smith n = mglevels[0]->levels; 1513aeaf226SBarry Smith for (i=0; i<n-1; i++) { 1523aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 1533aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 1543aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 1553aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 1563aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 1570de8493bSLawrence Mitchell ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr); 15873250ac0SBarry Smith ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 1593aeaf226SBarry Smith } 160391689abSStefano Zampini /* this is not null only if the smoother on the finest level 161391689abSStefano Zampini changes the rhs during PreSolve */ 162391689abSStefano Zampini ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr); 1633aeaf226SBarry Smith 1643aeaf226SBarry Smith for (i=0; i<n; i++) { 165*f3b08a26SMatthew G. Knepley if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {ierr = VecDestroy(&mglevels[i]->coarseSpace[c]);CHKERRQ(ierr);} 166*f3b08a26SMatthew G. Knepley ierr = PetscFree(mglevels[i]->coarseSpace);CHKERRQ(ierr); 167*f3b08a26SMatthew G. Knepley mglevels[i]->coarseSpace = NULL; 1683aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 1693aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1703aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 1713aeaf226SBarry Smith } 1723aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 1733aeaf226SBarry Smith } 174*f3b08a26SMatthew G. Knepley mg->Nc = 0; 1753aeaf226SBarry Smith } 1763aeaf226SBarry Smith PetscFunctionReturn(0); 1773aeaf226SBarry Smith } 1783aeaf226SBarry Smith 17997d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms) 1804b9ad928SBarry Smith { 181dfbe8321SBarry Smith PetscErrorCode ierr; 182f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 183ce94432eSBarry Smith MPI_Comm comm; 1843aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 18510eca3edSLisandro Dalcin PCMGType mgtype = mg->am; 18610167fecSBarry Smith PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V; 187f3fbd535SBarry Smith PetscInt i; 188f3fbd535SBarry Smith PetscMPIInt size; 189f3fbd535SBarry Smith const char *prefix; 190f3fbd535SBarry Smith PC ipc; 1913aeaf226SBarry Smith PetscInt n; 1924b9ad928SBarry Smith 1934b9ad928SBarry Smith PetscFunctionBegin; 1940700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 195548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 196548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 1971a97d7b6SStefano Zampini ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 1983aeaf226SBarry Smith if (mglevels) { 19910eca3edSLisandro Dalcin mgctype = mglevels[0]->cycles; 2003aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 2013aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 2023aeaf226SBarry Smith n = mglevels[0]->levels; 2033aeaf226SBarry Smith for (i=0; i<n; i++) { 2043aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2053aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 2063aeaf226SBarry Smith } 2073aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 2083aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 2093aeaf226SBarry Smith } 2103aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 2113aeaf226SBarry Smith } 212f3fbd535SBarry Smith 213f3fbd535SBarry Smith mg->nlevels = levels; 214f3fbd535SBarry Smith 215785e854fSJed Brown ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr); 2163bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 217f3fbd535SBarry Smith 218f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 219f3fbd535SBarry Smith 220a9db87a2SMatthew G Knepley mg->stageApply = 0; 221f3fbd535SBarry Smith for (i=0; i<levels; i++) { 222b00a9115SJed Brown ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr); 2232fa5cd67SKarl Rupp 224f3fbd535SBarry Smith mglevels[i]->level = i; 225f3fbd535SBarry Smith mglevels[i]->levels = levels; 22610eca3edSLisandro Dalcin mglevels[i]->cycles = mgctype; 227336babb1SJed Brown mg->default_smoothu = 2; 228336babb1SJed Brown mg->default_smoothd = 2; 22963e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 23063e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 23163e6d426SJed Brown mglevels[i]->eventresidual = 0; 23263e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 233f3fbd535SBarry Smith 234f3fbd535SBarry Smith if (comms) comm = comms[i]; 235f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 236422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr); 2370ee9a668SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 2380ee9a668SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 2390ee9a668SBarry Smith ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 2400ee9a668SBarry Smith if (i || levels == 1) { 2410ee9a668SBarry Smith char tprefix[128]; 2420ee9a668SBarry Smith 243336babb1SJed Brown ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 2440059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 245336babb1SJed Brown ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 246336babb1SJed Brown ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 247336babb1SJed Brown ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 2480ee9a668SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 249f3fbd535SBarry Smith 2500ee9a668SBarry Smith sprintf(tprefix,"mg_levels_%d_",(int)i); 2510ee9a668SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 2520ee9a668SBarry Smith } else { 253f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 254f3fbd535SBarry Smith 2559dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 256f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 257f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 258f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 259f3fbd535SBarry Smith if (size > 1) { 260f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 261f3fbd535SBarry Smith } else { 262f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 263f3fbd535SBarry Smith } 264753b7fb9SBarry Smith ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 265f3fbd535SBarry Smith } 2663bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr); 2672fa5cd67SKarl Rupp 268f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 26931567311SBarry Smith mg->rtol = 0.0; 27031567311SBarry Smith mg->abstol = 0.0; 27131567311SBarry Smith mg->dtol = 0.0; 27231567311SBarry Smith mg->ttol = 0.0; 27331567311SBarry Smith mg->cyclesperpcapply = 1; 274f3fbd535SBarry Smith } 275f3fbd535SBarry Smith mg->levels = mglevels; 27610eca3edSLisandro Dalcin ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 2774b9ad928SBarry Smith PetscFunctionReturn(0); 2784b9ad928SBarry Smith } 2794b9ad928SBarry Smith 28097d33e41SMatthew G. Knepley /*@C 28197d33e41SMatthew G. Knepley PCMGSetLevels - Sets the number of levels to use with MG. 28297d33e41SMatthew G. Knepley Must be called before any other MG routine. 28397d33e41SMatthew G. Knepley 28497d33e41SMatthew G. Knepley Logically Collective on PC 28597d33e41SMatthew G. Knepley 28697d33e41SMatthew G. Knepley Input Parameters: 28797d33e41SMatthew G. Knepley + pc - the preconditioner context 28897d33e41SMatthew G. Knepley . levels - the number of levels 28997d33e41SMatthew G. Knepley - comms - optional communicators for each level; this is to allow solving the coarser problems 29097d33e41SMatthew G. Knepley on smaller sets of processors. 29197d33e41SMatthew G. Knepley 29297d33e41SMatthew G. Knepley Level: intermediate 29397d33e41SMatthew G. Knepley 29497d33e41SMatthew G. Knepley Notes: 29597d33e41SMatthew G. Knepley If the number of levels is one then the multigrid uses the -mg_levels prefix 29697d33e41SMatthew G. Knepley for setting the level options rather than the -mg_coarse prefix. 29797d33e41SMatthew G. Knepley 29897d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels() 29997d33e41SMatthew G. Knepley @*/ 30097d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 30197d33e41SMatthew G. Knepley { 30297d33e41SMatthew G. Knepley PetscErrorCode ierr; 30397d33e41SMatthew G. Knepley 30497d33e41SMatthew G. Knepley PetscFunctionBegin; 30597d33e41SMatthew G. Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 30697d33e41SMatthew G. Knepley if (comms) PetscValidPointer(comms,3); 30737b5128cSMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr); 30897d33e41SMatthew G. Knepley PetscFunctionReturn(0); 30997d33e41SMatthew G. Knepley } 31097d33e41SMatthew G. Knepley 311c07bf074SBarry Smith 312c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 313c07bf074SBarry Smith { 314c07bf074SBarry Smith PetscErrorCode ierr; 315a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 316a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 317a06653b4SBarry Smith PetscInt i,n; 318c07bf074SBarry Smith 319c07bf074SBarry Smith PetscFunctionBegin; 320a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 321a06653b4SBarry Smith if (mglevels) { 322a06653b4SBarry Smith n = mglevels[0]->levels; 323a06653b4SBarry Smith for (i=0; i<n; i++) { 324a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 3256bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 326a06653b4SBarry Smith } 3276bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 328a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 329a06653b4SBarry Smith } 330a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 331a06653b4SBarry Smith } 332c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 333fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);CHKERRQ(ierr); 334fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);CHKERRQ(ierr); 335f3fbd535SBarry Smith PetscFunctionReturn(0); 336f3fbd535SBarry Smith } 337f3fbd535SBarry Smith 338f3fbd535SBarry Smith 339f3fbd535SBarry Smith 34009573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 34109573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 34209573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 343f3fbd535SBarry Smith 344f3fbd535SBarry Smith /* 345f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 346f3fbd535SBarry Smith or full cycle of multigrid. 347f3fbd535SBarry Smith 348f3fbd535SBarry Smith Note: 349f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 350f3fbd535SBarry Smith */ 351f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 352f3fbd535SBarry Smith { 353f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 354f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 355f3fbd535SBarry Smith PetscErrorCode ierr; 356391689abSStefano Zampini PC tpc; 357f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 358391689abSStefano Zampini PetscBool changeu,changed; 359f3fbd535SBarry Smith 360f3fbd535SBarry Smith PetscFunctionBegin; 361a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 362e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 363298cc208SBarry Smith for (i=0; i<levels; i++) { 364e1d8e5deSBarry Smith if (!mglevels[i]->A) { 36523ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 366298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 367e1d8e5deSBarry Smith } 368e1d8e5deSBarry Smith } 369e1d8e5deSBarry Smith 370391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 371391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 372391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 373391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 374391689abSStefano Zampini if (!changeu && !changed) { 375391689abSStefano Zampini ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 376f3fbd535SBarry Smith mglevels[levels-1]->b = b; 377391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 378391689abSStefano Zampini if (!mglevels[levels-1]->b) { 379391689abSStefano Zampini Vec *vec; 380391689abSStefano Zampini 381391689abSStefano Zampini ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 382391689abSStefano Zampini mglevels[levels-1]->b = *vec; 3837ae21d43SStefano Zampini ierr = PetscFree(vec);CHKERRQ(ierr); 384391689abSStefano Zampini } 385391689abSStefano Zampini ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 386391689abSStefano Zampini } 387f3fbd535SBarry Smith mglevels[levels-1]->x = x; 388391689abSStefano Zampini 38931567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 390f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 39131567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 3920298fd71SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr); 393f3fbd535SBarry Smith } 3942fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 39531567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 3962fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 39731567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 3982fa5cd67SKarl Rupp } else { 399f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 400f3fbd535SBarry Smith } 401a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 402391689abSStefano Zampini if (!changeu && !changed) mglevels[levels-1]->b = NULL; 403f3fbd535SBarry Smith PetscFunctionReturn(0); 404f3fbd535SBarry Smith } 405f3fbd535SBarry Smith 406f3fbd535SBarry Smith 4074416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 408f3fbd535SBarry Smith { 409f3fbd535SBarry Smith PetscErrorCode ierr; 410710315b6SLawrence Mitchell PetscInt levels,cycles; 411*f3b08a26SMatthew G. Knepley PetscBool flg, flg2; 412f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 4133d3eaba7SBarry Smith PC_MG_Levels **mglevels; 414f3fbd535SBarry Smith PCMGType mgtype; 415f3fbd535SBarry Smith PCMGCycleType mgctype; 4162134b1e4SBarry Smith PCMGGalerkinType gtype; 417f3fbd535SBarry Smith 418f3fbd535SBarry Smith PetscFunctionBegin; 4191d743356SStefano Zampini levels = PetscMax(mg->nlevels,1); 420e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr); 421f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 4221a97d7b6SStefano Zampini if (!flg && !mg->levels && pc->dm) { 423eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 424eb3f98d2SBarry Smith levels++; 4253aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 426eb3f98d2SBarry Smith } 4270298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 4283aeaf226SBarry Smith mglevels = mg->levels; 4293aeaf226SBarry Smith 430f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 431f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 432f3fbd535SBarry Smith if (flg) { 433f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 4342fa5cd67SKarl Rupp } 4352134b1e4SBarry Smith gtype = mg->galerkin; 4362134b1e4SBarry Smith ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg);CHKERRQ(ierr); 4372134b1e4SBarry Smith if (flg) { 4382134b1e4SBarry Smith ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr); 439f3fbd535SBarry Smith } 440*f3b08a26SMatthew G. Knepley flg2 = PETSC_FALSE; 441*f3b08a26SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr); 442*f3b08a26SMatthew G. Knepley if (flg) {ierr = PCMGSetAdaptInterpolation(pc, flg2);CHKERRQ(ierr);} 443*f3b08a26SMatthew G. Knepley ierr = PetscOptionsInt("-pc_mg_adapt_interp_n","Size of the coarse space for adaptive interpolation","PCMGSetCoarseSpace",mg->Nc,&mg->Nc,&flg);CHKERRQ(ierr); 444*f3b08a26SMatthew G. Knepley ierr = PetscOptionsEnum("-pc_mg_adapt_interp_coarse_space","Type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector","PCMGSetAdaptCoarseSpaceType",PCMGCoarseSpaceTypes,(PetscEnum)mg->coarseSpaceType,(PetscEnum*)&mg->coarseSpaceType,&flg);CHKERRQ(ierr); 445*f3b08a26SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);CHKERRQ(ierr); 446f56b1265SBarry Smith flg = PETSC_FALSE; 4478e5aa403SBarry Smith ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr); 448f442ab6aSBarry Smith if (flg) { 449f442ab6aSBarry Smith ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr); 450f442ab6aSBarry Smith } 45131567311SBarry Smith mgtype = mg->am; 452f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 453f3fbd535SBarry Smith if (flg) { 454f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 455f3fbd535SBarry Smith } 45631567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 4575363c1e0SLisandro Dalcin ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 458f3fbd535SBarry Smith if (flg) { 459f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 460f3fbd535SBarry Smith } 461f3fbd535SBarry Smith } 462f3fbd535SBarry Smith flg = PETSC_FALSE; 4630298fd71SBarry Smith ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 464f3fbd535SBarry Smith if (flg) { 465f3fbd535SBarry Smith PetscInt i; 466f3fbd535SBarry Smith char eventname[128]; 4671a97d7b6SStefano Zampini 468f3fbd535SBarry Smith levels = mglevels[0]->levels; 469f3fbd535SBarry Smith for (i=0; i<levels; i++) { 470f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 47163e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 472f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 47363e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 474f3fbd535SBarry Smith if (i) { 475f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 47663e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 477f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 47863e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 479f3fbd535SBarry Smith } 480f3fbd535SBarry Smith } 481eec35531SMatthew G Knepley 482ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 483eec35531SMatthew G Knepley { 484eec35531SMatthew G Knepley const char *sname = "MG Apply"; 485eec35531SMatthew G Knepley PetscStageLog stageLog; 486eec35531SMatthew G Knepley PetscInt st; 487eec35531SMatthew G Knepley 488eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 489eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 490eec35531SMatthew G Knepley PetscBool same; 491eec35531SMatthew G Knepley 492eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 4932fa5cd67SKarl Rupp if (same) mg->stageApply = st; 494eec35531SMatthew G Knepley } 495eec35531SMatthew G Knepley if (!mg->stageApply) { 496eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 497eec35531SMatthew G Knepley } 498eec35531SMatthew G Knepley } 499ec60ca73SSatish Balay #endif 500f3fbd535SBarry Smith } 501f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 502*f3b08a26SMatthew G. Knepley /* Check option consistency */ 503*f3b08a26SMatthew G. Knepley ierr = PCMGGetGalerkin(pc, >ype);CHKERRQ(ierr); 504*f3b08a26SMatthew G. Knepley ierr = PCMGGetAdaptInterpolation(pc, &flg);CHKERRQ(ierr); 505*f3b08a26SMatthew G. Knepley if (flg && (gtype >= PC_MG_GALERKIN_NONE)) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_INCOMP, "Must use Galerkin coarse operators when adapting the interpolator"); 506f3fbd535SBarry Smith PetscFunctionReturn(0); 507f3fbd535SBarry Smith } 508f3fbd535SBarry Smith 5090a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL}; 5100a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL}; 5110a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL}; 512*f3b08a26SMatthew G. Knepley const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL}; 513f3fbd535SBarry Smith 5149804daf3SBarry Smith #include <petscdraw.h> 515f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 516f3fbd535SBarry Smith { 517f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 518f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 519f3fbd535SBarry Smith PetscErrorCode ierr; 520e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 521219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 522f3fbd535SBarry Smith 523f3fbd535SBarry Smith PetscFunctionBegin; 524251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 5255b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 526219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 527f3fbd535SBarry Smith if (iascii) { 528e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 529efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 53031567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 53131567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 532f3fbd535SBarry Smith } 5332134b1e4SBarry Smith if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 534f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 5352134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 5362134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr); 5372134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 5382134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr); 5392134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 5402134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr); 5414f66f45eSBarry Smith } else { 5424f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 543f3fbd535SBarry Smith } 5445adeb434SBarry Smith if (mg->view){ 5455adeb434SBarry Smith ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr); 5465adeb434SBarry Smith } 547f3fbd535SBarry Smith for (i=0; i<levels; i++) { 548f3fbd535SBarry Smith if (!i) { 54989cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 550f3fbd535SBarry Smith } else { 55189cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 552f3fbd535SBarry Smith } 553f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 554f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 555f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 556f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 557f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 558f3fbd535SBarry Smith } else if (i) { 55989cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 560f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 561f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 562f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 563f3fbd535SBarry Smith } 564f3fbd535SBarry Smith } 5655b0b0462SBarry Smith } else if (isbinary) { 5665b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 5675b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 5685b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 5695b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 5705b0b0462SBarry Smith } 5715b0b0462SBarry Smith } 572219b1045SBarry Smith } else if (isdraw) { 573219b1045SBarry Smith PetscDraw draw; 574b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 575219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 576219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 5770298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 578219b1045SBarry Smith bottom = y - th; 579219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 580b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 581219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 582219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 583219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 584b4375e8dSPeter Brune } else { 585b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 586b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 587b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 588b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 589b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 590b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 591b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 592b4375e8dSPeter Brune } 5930298fd71SBarry Smith ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 5941cd381d6SBarry Smith bottom -= th; 595219b1045SBarry Smith } 5965b0b0462SBarry Smith } 597f3fbd535SBarry Smith PetscFunctionReturn(0); 598f3fbd535SBarry Smith } 599f3fbd535SBarry Smith 600af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 601af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 602cab2e9ccSBarry Smith 603f3fbd535SBarry Smith /* 604f3fbd535SBarry Smith Calls setup for the KSP on each level 605f3fbd535SBarry Smith */ 606f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 607f3fbd535SBarry Smith { 608f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 609f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 610f3fbd535SBarry Smith PetscErrorCode ierr; 6111a97d7b6SStefano Zampini PetscInt i,n; 61298e04984SBarry Smith PC cpc; 6133db492dfSBarry Smith PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 614f3fbd535SBarry Smith Mat dA,dB; 615f3fbd535SBarry Smith Vec tvec; 616218a07d4SBarry Smith DM *dms; 6170a545947SLisandro Dalcin PetscViewer viewer = NULL; 6182134b1e4SBarry Smith PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE; 619f3fbd535SBarry Smith 620f3fbd535SBarry Smith PetscFunctionBegin; 6211a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up"); 6221a97d7b6SStefano Zampini n = mglevels[0]->levels; 62301bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 6243aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 6253aeaf226SBarry Smith PetscInt levels; 6263aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 6273aeaf226SBarry Smith levels++; 6283aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 6290298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 6303aeaf226SBarry Smith n = levels; 6313aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 6323aeaf226SBarry Smith mglevels = mg->levels; 6333aeaf226SBarry Smith } 6343aeaf226SBarry Smith } 63598e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 6363aeaf226SBarry Smith 637f3fbd535SBarry Smith 638f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 639f3fbd535SBarry Smith /* so use those from global PC */ 640f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 6410298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 64298e04984SBarry Smith if (opsset) { 64398e04984SBarry Smith Mat mmat; 64423ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 64598e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 64698e04984SBarry Smith } 6474a5f07a5SMark F. Adams 6484a5f07a5SMark F. Adams if (!opsset) { 64971b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 6504a5f07a5SMark F. Adams if (use_amat) { 651f3fbd535SBarry 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); 65223ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 65369858f1bSStefano Zampini } else { 6544a5f07a5SMark 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); 65523ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 6564a5f07a5SMark F. Adams } 6574a5f07a5SMark F. Adams } 658f3fbd535SBarry Smith 6599c7ffb39SBarry Smith for (i=n-1; i>0; i--) { 6609c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 6619c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 6629c7ffb39SBarry Smith continue; 6639c7ffb39SBarry Smith } 6649c7ffb39SBarry Smith } 6652134b1e4SBarry Smith 6662134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 6672134b1e4SBarry Smith if (dA == dB) dAeqdB = PETSC_TRUE; 6682134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) { 6692134b1e4SBarry Smith needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 6702134b1e4SBarry Smith } 6712134b1e4SBarry Smith 6722134b1e4SBarry Smith 6739c7ffb39SBarry Smith /* 6749c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 6759c7ffb39SBarry Smith Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 6769c7ffb39SBarry Smith */ 6772134b1e4SBarry Smith if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 6782d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 6792e499ae9SBarry Smith Mat p; 68073250ac0SBarry Smith Vec rscale; 681785e854fSJed Brown ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 682218a07d4SBarry Smith dms[n-1] = pc->dm; 683ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 684ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 68541fce8fdSHong Zhang /* 68641fce8fdSHong Zhang Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level. 68741fce8fdSHong Zhang Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options. 68841fce8fdSHong Zhang But it is safe to use -dm_mat_type. 68941fce8fdSHong Zhang 69041fce8fdSHong Zhang The mat type should not be hardcoded like this, we need to find a better way. 69141fce8fdSHong Zhang ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr); 69241fce8fdSHong Zhang */ 693218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 694942e3340SBarry Smith DMKSP kdm; 695eab52d2bSLawrence Mitchell PetscBool dmhasrestrict, dmhasinject; 6963c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 6972134b1e4SBarry Smith if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 698c27ee7a3SPatrick Farrell if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 699c27ee7a3SPatrick Farrell ierr = KSPSetDM(mglevels[i]->smoothu,dms[i]);CHKERRQ(ierr); 700c27ee7a3SPatrick Farrell if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);CHKERRQ(ierr);} 701c27ee7a3SPatrick Farrell } 702942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 703d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 704d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 7050298fd71SBarry Smith kdm->ops->computerhs = NULL; 7060298fd71SBarry Smith kdm->rhsctx = NULL; 70724c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 708e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 7092d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 71000ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 71173250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 7126bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 713218a07d4SBarry Smith } 7143ad4599aSBarry Smith ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr); 7153ad4599aSBarry Smith if (dmhasrestrict && !mglevels[i+1]->restrct){ 7163ad4599aSBarry Smith ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr); 7173ad4599aSBarry Smith ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr); 7183ad4599aSBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 7193ad4599aSBarry Smith } 720eab52d2bSLawrence Mitchell ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr); 721eab52d2bSLawrence Mitchell if (dmhasinject && !mglevels[i+1]->inject){ 722eab52d2bSLawrence Mitchell ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr); 723eab52d2bSLawrence Mitchell ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr); 724eab52d2bSLawrence Mitchell ierr = MatDestroy(&p);CHKERRQ(ierr); 725eab52d2bSLawrence Mitchell } 72624c3aa18SBarry Smith } 7272d2b81a6SBarry Smith 728ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 7292d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 730ce4cda84SJed Brown } 731cab2e9ccSBarry Smith 732ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 7332134b1e4SBarry Smith /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 734cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 735cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 736c27ee7a3SPatrick Farrell if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) { 737c27ee7a3SPatrick Farrell ierr = KSPSetDM(mglevels[n-1]->smoothu,pc->dm);CHKERRQ(ierr); 738c27ee7a3SPatrick Farrell ierr = KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);CHKERRQ(ierr); 739c27ee7a3SPatrick Farrell } 740218a07d4SBarry Smith } 741218a07d4SBarry Smith 7422134b1e4SBarry Smith if (mg->galerkin < PC_MG_GALERKIN_NONE) { 7432134b1e4SBarry Smith Mat A,B; 7442134b1e4SBarry Smith PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 7452134b1e4SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 7462134b1e4SBarry Smith 7472134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 7482134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE; 7492134b1e4SBarry Smith if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 750f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 751b9367914SBarry 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"); 752b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 753b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 754b9367914SBarry Smith } 755b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 756b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 757b9367914SBarry Smith } 7582134b1e4SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 7592134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr); 7602134b1e4SBarry Smith } 7612134b1e4SBarry Smith if (doA) { 7622df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr); 7632134b1e4SBarry Smith } 7642134b1e4SBarry Smith if (doB) { 7652df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr); 766b9367914SBarry Smith } 7672134b1e4SBarry Smith /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 7682134b1e4SBarry Smith if (!doA && dAeqdB) { 7692134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);} 7702134b1e4SBarry Smith A = B; 7712134b1e4SBarry Smith } else if (!doA && reuse == MAT_INITIAL_MATRIX ) { 7722134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr); 7732134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 774b9367914SBarry Smith } 7752134b1e4SBarry Smith if (!doB && dAeqdB) { 7762134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);} 7772134b1e4SBarry Smith B = A; 7782134b1e4SBarry Smith } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 77923ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 7802134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 7817d537d92SJed Brown } 7822134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 7832134b1e4SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr); 7842134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr); 7852134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr); 7862134b1e4SBarry Smith } 7872134b1e4SBarry Smith dA = A; 788f3fbd535SBarry Smith dB = B; 789f3fbd535SBarry Smith } 790f3fbd535SBarry Smith } 791*f3b08a26SMatthew G. Knepley 792*f3b08a26SMatthew G. Knepley /* Create key for storing eigenvalues inside eigenvectors */ 793*f3b08a26SMatthew G. Knepley if (mg->adaptInterpolation && mg->coarseSpaceType == PCMG_EIGENVECTOR && mg->eigenvalue < 0) {ierr = PetscObjectComposedDataRegister(&mg->eigenvalue);CHKERRQ(ierr);} 794*f3b08a26SMatthew G. Knepley 795*f3b08a26SMatthew G. Knepley /* Adapt interpolation matrices */ 796*f3b08a26SMatthew G. Knepley if (mg->adaptInterpolation) { 797*f3b08a26SMatthew G. Knepley mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */ 798*f3b08a26SMatthew G. Knepley for (i = 0; i < n; ++i) { 799*f3b08a26SMatthew G. Knepley ierr = PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);CHKERRQ(ierr); 800*f3b08a26SMatthew G. Knepley if (i) {ierr = PCMGAdaptInterpolator_Internal(pc, i, mglevels[i-1]->smoothu, mglevels[i]->smoothu, mg->Nc, mglevels[i-1]->coarseSpace, mglevels[i]->coarseSpace);CHKERRQ(ierr);} 801*f3b08a26SMatthew G. Knepley } 802*f3b08a26SMatthew G. Knepley for (i = n-2; i > -1; --i) { 803*f3b08a26SMatthew G. Knepley ierr = PCMGRecomputeLevelOperators_Internal(pc, i);CHKERRQ(ierr); 804*f3b08a26SMatthew G. Knepley } 805*f3b08a26SMatthew G. Knepley } 806*f3b08a26SMatthew G. Knepley 8072134b1e4SBarry Smith if (needRestricts && pc->dm) { 808caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 809ccceb50cSJed Brown DM dmfine,dmcoarse; 810caa4e7f2SJed Brown Mat Restrict,Inject; 811caa4e7f2SJed Brown Vec rscale; 812ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 813ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 814caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 815caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 816eab52d2bSLawrence Mitchell ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr); 817caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 818caa4e7f2SJed Brown } 819caa4e7f2SJed Brown } 820f3fbd535SBarry Smith 821f3fbd535SBarry Smith if (!pc->setupcalled) { 822f3fbd535SBarry Smith for (i=0; i<n; i++) { 823f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 824f3fbd535SBarry Smith } 825f3fbd535SBarry Smith for (i=1; i<n; i++) { 826f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 827f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 828f3fbd535SBarry Smith } 829f3fbd535SBarry Smith } 8303ad4599aSBarry Smith /* insure that if either interpolation or restriction is set the other other one is set */ 831f3fbd535SBarry Smith for (i=1; i<n; i++) { 8323ad4599aSBarry Smith ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr); 8333ad4599aSBarry Smith ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr); 834f3fbd535SBarry Smith } 835f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 836f3fbd535SBarry Smith if (!mglevels[i]->b) { 837f3fbd535SBarry Smith Vec *vec; 8382a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 839f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 8406bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 841f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 842f3fbd535SBarry Smith } 843f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 844f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 845f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 8466bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 847f3fbd535SBarry Smith } 848f3fbd535SBarry Smith if (!mglevels[i]->x) { 849f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 850f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 8516bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 852f3fbd535SBarry Smith } 853f3fbd535SBarry Smith } 854f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 855f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 856f3fbd535SBarry Smith Vec *vec; 8572a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 858f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 8596bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 860f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 861f3fbd535SBarry Smith } 862f3fbd535SBarry Smith } 863f3fbd535SBarry Smith 86498e04984SBarry Smith if (pc->dm) { 86598e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 86698e04984SBarry Smith for (i=0; i<n-1; i++) { 86798e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 86898e04984SBarry Smith } 86998e04984SBarry Smith } 870f3fbd535SBarry Smith 871f3fbd535SBarry Smith for (i=1; i<n; i++) { 8722a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 873f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 874f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 875f3fbd535SBarry Smith } 87663e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 877f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 878c0decd05SBarry Smith if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 879899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 880899639b0SHong Zhang } 88163e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 882d42688cbSBarry Smith if (!mglevels[i]->residual) { 883d42688cbSBarry Smith Mat mat; 88413842ffbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 88554b2cd4bSJed Brown ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 886d42688cbSBarry Smith } 887f3fbd535SBarry Smith } 888f3fbd535SBarry Smith for (i=1; i<n; i++) { 889f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 890f3fbd535SBarry Smith Mat downmat,downpmat; 891f3fbd535SBarry Smith 892f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 8930298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 894f3fbd535SBarry Smith if (!opsset) { 89523ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 89623ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 897f3fbd535SBarry Smith } 898f3fbd535SBarry Smith 899f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 90063e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 901f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 902c0decd05SBarry Smith if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) { 903899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 904899639b0SHong Zhang } 90563e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 906f3fbd535SBarry Smith } 907f3fbd535SBarry Smith } 908f3fbd535SBarry Smith 90963e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 910f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 911c0decd05SBarry Smith if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 912899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 913899639b0SHong Zhang } 91463e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 915f3fbd535SBarry Smith 916*f3b08a26SMatthew G. Knepley /* Adapt interpolation matrices */ 917*f3b08a26SMatthew G. Knepley if (mg->adaptInterpolation) { 918*f3b08a26SMatthew G. Knepley Mat InterpAdapt; 919*f3b08a26SMatthew G. Knepley 920*f3b08a26SMatthew G. Knepley /* Create coarse modes */ 921*f3b08a26SMatthew G. Knepley if (!mglevels[0]->coarseSpace) { 922*f3b08a26SMatthew G. Knepley DM dm; 923*f3b08a26SMatthew G. Knepley mg->Nc = mg->Nc ? mg->Nc : 6; 924*f3b08a26SMatthew G. Knepley for (i = 0; i < n; ++i) { 925*f3b08a26SMatthew G. Knepley ierr = KSPGetDM(mglevels[i]->smoothu, &dm);CHKERRQ(ierr); 926*f3b08a26SMatthew G. Knepley ierr = PCMGCreateCoarseSpace_Private(pc, dm, mglevels[i]->smoothu, mg->Nc, &mglevels[i]->coarseSpace);CHKERRQ(ierr); 927*f3b08a26SMatthew G. Knepley } 928*f3b08a26SMatthew G. Knepley } 929*f3b08a26SMatthew G. Knepley /* Adapt interpolators */ 930*f3b08a26SMatthew G. Knepley for (i = 1; i < n; ++i) { 931*f3b08a26SMatthew G. Knepley DM dm, cdm; 932*f3b08a26SMatthew G. Knepley 933*f3b08a26SMatthew G. Knepley ierr = KSPGetDM(mglevels[i-1]->smoothu, &cdm);CHKERRQ(ierr); 934*f3b08a26SMatthew G. Knepley ierr = KSPGetDM(mglevels[i]->smoothu, &dm);CHKERRQ(ierr); 935*f3b08a26SMatthew G. Knepley ierr = DMAdaptInterpolator(cdm, dm, mglevels[i]->interpolate, mglevels[i]->smoothu, mg->Nc, mglevels[i]->coarseSpace, mglevels[i-1]->coarseSpace, &InterpAdapt, pc);CHKERRQ(ierr); 936*f3b08a26SMatthew G. Knepley ierr = PCMGSetInterpolation(pc, i, InterpAdapt);CHKERRQ(ierr); 937*f3b08a26SMatthew G. Knepley ierr = MatDestroy(&InterpAdapt);CHKERRQ(ierr); 938*f3b08a26SMatthew G. Knepley } 939*f3b08a26SMatthew G. Knepley } 940*f3b08a26SMatthew G. Knepley 941f3fbd535SBarry Smith /* 942f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 943e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 944f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 945f3fbd535SBarry Smith 946f3fbd535SBarry Smith Only support one or the other at the same time. 947f3fbd535SBarry Smith */ 948f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 949c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 950ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 951f3fbd535SBarry Smith dump = PETSC_FALSE; 952f3fbd535SBarry Smith #endif 953c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 954ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 955f3fbd535SBarry Smith 956f3fbd535SBarry Smith if (viewer) { 957f3fbd535SBarry Smith for (i=1; i<n; i++) { 958f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 959f3fbd535SBarry Smith } 960f3fbd535SBarry Smith for (i=0; i<n; i++) { 961f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 962f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 963f3fbd535SBarry Smith } 964f3fbd535SBarry Smith } 965f3fbd535SBarry Smith PetscFunctionReturn(0); 966f3fbd535SBarry Smith } 967f3fbd535SBarry Smith 968f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 969f3fbd535SBarry Smith 97097d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 97197d33e41SMatthew G. Knepley { 97297d33e41SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 97397d33e41SMatthew G. Knepley 97497d33e41SMatthew G. Knepley PetscFunctionBegin; 97597d33e41SMatthew G. Knepley *levels = mg->nlevels; 97697d33e41SMatthew G. Knepley PetscFunctionReturn(0); 97797d33e41SMatthew G. Knepley } 97897d33e41SMatthew G. Knepley 9794b9ad928SBarry Smith /*@ 98097177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 9814b9ad928SBarry Smith 9824b9ad928SBarry Smith Not Collective 9834b9ad928SBarry Smith 9844b9ad928SBarry Smith Input Parameter: 9854b9ad928SBarry Smith . pc - the preconditioner context 9864b9ad928SBarry Smith 9874b9ad928SBarry Smith Output parameter: 9884b9ad928SBarry Smith . levels - the number of levels 9894b9ad928SBarry Smith 9904b9ad928SBarry Smith Level: advanced 9914b9ad928SBarry Smith 99297177400SBarry Smith .seealso: PCMGSetLevels() 9934b9ad928SBarry Smith @*/ 9947087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 9954b9ad928SBarry Smith { 996e952c529SMatthew G. Knepley PetscErrorCode ierr; 9974b9ad928SBarry Smith 9984b9ad928SBarry Smith PetscFunctionBegin; 9990700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10004482741eSBarry Smith PetscValidIntPointer(levels,2); 100197d33e41SMatthew G. Knepley *levels = 0; 100237b5128cSMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr); 10034b9ad928SBarry Smith PetscFunctionReturn(0); 10044b9ad928SBarry Smith } 10054b9ad928SBarry Smith 10064b9ad928SBarry Smith /*@ 100797177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 10084b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 10094b9ad928SBarry Smith 1010ad4df100SBarry Smith Logically Collective on PC 10114b9ad928SBarry Smith 10124b9ad928SBarry Smith Input Parameters: 10134b9ad928SBarry Smith + pc - the preconditioner context 10149dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 10159dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 10164b9ad928SBarry Smith 10174b9ad928SBarry Smith Options Database Key: 10184b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 10194b9ad928SBarry Smith additive, full, kaskade 10204b9ad928SBarry Smith 10214b9ad928SBarry Smith Level: advanced 10224b9ad928SBarry Smith 102397177400SBarry Smith .seealso: PCMGSetLevels() 10244b9ad928SBarry Smith @*/ 10257087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 10264b9ad928SBarry Smith { 1027f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 10284b9ad928SBarry Smith 10294b9ad928SBarry Smith PetscFunctionBegin; 10300700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1031c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 103231567311SBarry Smith mg->am = form; 10339dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 1034c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 1035c60c7ad4SBarry Smith PetscFunctionReturn(0); 1036c60c7ad4SBarry Smith } 1037c60c7ad4SBarry Smith 1038c60c7ad4SBarry Smith /*@ 1039c60c7ad4SBarry Smith PCMGGetType - Determines the form of multigrid to use: 1040c60c7ad4SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 1041c60c7ad4SBarry Smith 1042c60c7ad4SBarry Smith Logically Collective on PC 1043c60c7ad4SBarry Smith 1044c60c7ad4SBarry Smith Input Parameter: 1045c60c7ad4SBarry Smith . pc - the preconditioner context 1046c60c7ad4SBarry Smith 1047c60c7ad4SBarry Smith Output Parameter: 1048c60c7ad4SBarry Smith . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 1049c60c7ad4SBarry Smith 1050c60c7ad4SBarry Smith 1051c60c7ad4SBarry Smith Level: advanced 1052c60c7ad4SBarry Smith 1053c60c7ad4SBarry Smith .seealso: PCMGSetLevels() 1054c60c7ad4SBarry Smith @*/ 1055c60c7ad4SBarry Smith PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 1056c60c7ad4SBarry Smith { 1057c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1058c60c7ad4SBarry Smith 1059c60c7ad4SBarry Smith PetscFunctionBegin; 1060c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1061c60c7ad4SBarry Smith *type = mg->am; 10624b9ad928SBarry Smith PetscFunctionReturn(0); 10634b9ad928SBarry Smith } 10644b9ad928SBarry Smith 10654b9ad928SBarry Smith /*@ 10660d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 10674b9ad928SBarry Smith complicated cycling. 10684b9ad928SBarry Smith 1069ad4df100SBarry Smith Logically Collective on PC 10704b9ad928SBarry Smith 10714b9ad928SBarry Smith Input Parameters: 1072c2be2410SBarry Smith + pc - the multigrid context 1073c1cbb1deSBarry Smith - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 10744b9ad928SBarry Smith 10754b9ad928SBarry Smith Options Database Key: 1076c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired 10774b9ad928SBarry Smith 10784b9ad928SBarry Smith Level: advanced 10794b9ad928SBarry Smith 10800d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 10814b9ad928SBarry Smith @*/ 10827087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 10834b9ad928SBarry Smith { 1084f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1085f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 108679416396SBarry Smith PetscInt i,levels; 10874b9ad928SBarry Smith 10884b9ad928SBarry Smith PetscFunctionBegin; 10890700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 109069fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc,n,2); 10911a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1092f3fbd535SBarry Smith levels = mglevels[0]->levels; 10932fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 10944b9ad928SBarry Smith PetscFunctionReturn(0); 10954b9ad928SBarry Smith } 10964b9ad928SBarry Smith 10978cc2d5dfSBarry Smith /*@ 10988cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 10998cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 11008cc2d5dfSBarry Smith 1101ad4df100SBarry Smith Logically Collective on PC 11028cc2d5dfSBarry Smith 11038cc2d5dfSBarry Smith Input Parameters: 11048cc2d5dfSBarry Smith + pc - the multigrid context 11058cc2d5dfSBarry Smith - n - number of cycles (default is 1) 11068cc2d5dfSBarry Smith 11078cc2d5dfSBarry Smith Options Database Key: 1108e1bc860dSBarry Smith . -pc_mg_multiplicative_cycles n 11098cc2d5dfSBarry Smith 11108cc2d5dfSBarry Smith Level: advanced 11118cc2d5dfSBarry Smith 111295452b02SPatrick Sanan Notes: 111395452b02SPatrick Sanan This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 11148cc2d5dfSBarry Smith 11158cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 11168cc2d5dfSBarry Smith @*/ 11177087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 11188cc2d5dfSBarry Smith { 1119f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 11208cc2d5dfSBarry Smith 11218cc2d5dfSBarry Smith PetscFunctionBegin; 11220700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1123c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 11242a21e185SBarry Smith mg->cyclesperpcapply = n; 11258cc2d5dfSBarry Smith PetscFunctionReturn(0); 11268cc2d5dfSBarry Smith } 11278cc2d5dfSBarry Smith 11282134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1129fb15c04eSBarry Smith { 1130fb15c04eSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1131fb15c04eSBarry Smith 1132fb15c04eSBarry Smith PetscFunctionBegin; 11332134b1e4SBarry Smith mg->galerkin = use; 1134fb15c04eSBarry Smith PetscFunctionReturn(0); 1135fb15c04eSBarry Smith } 1136fb15c04eSBarry Smith 1137c2be2410SBarry Smith /*@ 113897177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 113982b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1140c2be2410SBarry Smith 1141ad4df100SBarry Smith Logically Collective on PC 1142c2be2410SBarry Smith 1143c2be2410SBarry Smith Input Parameters: 1144c91913d3SJed Brown + pc - the multigrid context 11452134b1e4SBarry Smith - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1146c2be2410SBarry Smith 1147c2be2410SBarry Smith Options Database Key: 11482134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> 1149c2be2410SBarry Smith 1150c2be2410SBarry Smith Level: intermediate 1151c2be2410SBarry Smith 115295452b02SPatrick Sanan Notes: 115395452b02SPatrick Sanan Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 11541c1aac46SBarry Smith use the PCMG construction of the coarser grids. 11551c1aac46SBarry Smith 11562134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType 11573fc8bf9cSBarry Smith 1158c2be2410SBarry Smith @*/ 11592134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1160c2be2410SBarry Smith { 1161fb15c04eSBarry Smith PetscErrorCode ierr; 1162c2be2410SBarry Smith 1163c2be2410SBarry Smith PetscFunctionBegin; 11640700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11652134b1e4SBarry Smith ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr); 1166c2be2410SBarry Smith PetscFunctionReturn(0); 1167c2be2410SBarry Smith } 1168c2be2410SBarry Smith 11693fc8bf9cSBarry Smith /*@ 11703fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 117182b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 11723fc8bf9cSBarry Smith 11733fc8bf9cSBarry Smith Not Collective 11743fc8bf9cSBarry Smith 11753fc8bf9cSBarry Smith Input Parameter: 11763fc8bf9cSBarry Smith . pc - the multigrid context 11773fc8bf9cSBarry Smith 11783fc8bf9cSBarry Smith Output Parameter: 11792134b1e4SBarry 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 11803fc8bf9cSBarry Smith 11813fc8bf9cSBarry Smith Level: intermediate 11823fc8bf9cSBarry Smith 11832134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType 11843fc8bf9cSBarry Smith 11853fc8bf9cSBarry Smith @*/ 11862134b1e4SBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 11873fc8bf9cSBarry Smith { 1188f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 11893fc8bf9cSBarry Smith 11903fc8bf9cSBarry Smith PetscFunctionBegin; 11910700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11922134b1e4SBarry Smith *galerkin = mg->galerkin; 11933fc8bf9cSBarry Smith PetscFunctionReturn(0); 11943fc8bf9cSBarry Smith } 11953fc8bf9cSBarry Smith 1196*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt) 1197*f3b08a26SMatthew G. Knepley { 1198*f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1199*f3b08a26SMatthew G. Knepley 1200*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1201*f3b08a26SMatthew G. Knepley mg->adaptInterpolation = adapt; 1202*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1203*f3b08a26SMatthew G. Knepley } 1204*f3b08a26SMatthew G. Knepley 1205*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt) 1206*f3b08a26SMatthew G. Knepley { 1207*f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1208*f3b08a26SMatthew G. Knepley 1209*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1210*f3b08a26SMatthew G. Knepley *adapt = mg->adaptInterpolation; 1211*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1212*f3b08a26SMatthew G. Knepley } 1213*f3b08a26SMatthew G. Knepley 1214*f3b08a26SMatthew G. Knepley /*@ 1215*f3b08a26SMatthew G. Knepley PCMGSetAdaptInterpolation - Adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1216*f3b08a26SMatthew G. Knepley 1217*f3b08a26SMatthew G. Knepley Logically Collective on PC 1218*f3b08a26SMatthew G. Knepley 1219*f3b08a26SMatthew G. Knepley Input Parameters: 1220*f3b08a26SMatthew G. Knepley + pc - the multigrid context 1221*f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator 1222*f3b08a26SMatthew G. Knepley 1223*f3b08a26SMatthew G. Knepley Options Database Keys: 1224*f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp - Turn on adaptation 1225*f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension 1226*f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector 1227*f3b08a26SMatthew G. Knepley 1228*f3b08a26SMatthew G. Knepley Level: intermediate 1229*f3b08a26SMatthew G. Knepley 1230*f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1231*f3b08a26SMatthew G. Knepley .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin() 1232*f3b08a26SMatthew G. Knepley @*/ 1233*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt) 1234*f3b08a26SMatthew G. Knepley { 1235*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1236*f3b08a26SMatthew G. Knepley 1237*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1238*f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1239*f3b08a26SMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));CHKERRQ(ierr); 1240*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1241*f3b08a26SMatthew G. Knepley } 1242*f3b08a26SMatthew G. Knepley 1243*f3b08a26SMatthew G. Knepley /*@ 1244*f3b08a26SMatthew G. Knepley PCMGGetAdaptInterpolation - Get the flag to adapt the interpolator based upon a vector space which should be accurately captured by the next coarser mesh, and thus accurately interpolated. 1245*f3b08a26SMatthew G. Knepley 1246*f3b08a26SMatthew G. Knepley Not collective 1247*f3b08a26SMatthew G. Knepley 1248*f3b08a26SMatthew G. Knepley Input Parameter: 1249*f3b08a26SMatthew G. Knepley . pc - the multigrid context 1250*f3b08a26SMatthew G. Knepley 1251*f3b08a26SMatthew G. Knepley Output Parameter: 1252*f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator 1253*f3b08a26SMatthew G. Knepley 1254*f3b08a26SMatthew G. Knepley Level: intermediate 1255*f3b08a26SMatthew G. Knepley 1256*f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1257*f3b08a26SMatthew G. Knepley .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin() 1258*f3b08a26SMatthew G. Knepley @*/ 1259*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt) 1260*f3b08a26SMatthew G. Knepley { 1261*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1262*f3b08a26SMatthew G. Knepley 1263*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1264*f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1265*f3b08a26SMatthew G. Knepley PetscValidPointer(adapt, 2); 1266*f3b08a26SMatthew G. Knepley ierr = PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));CHKERRQ(ierr); 1267*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1268*f3b08a26SMatthew G. Knepley } 1269*f3b08a26SMatthew G. Knepley 12704b9ad928SBarry Smith /*@ 127106792cafSBarry Smith PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1272710315b6SLawrence Mitchell on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1273710315b6SLawrence Mitchell pre- and post-smoothing steps. 127406792cafSBarry Smith 127506792cafSBarry Smith Logically Collective on PC 127606792cafSBarry Smith 127706792cafSBarry Smith Input Parameters: 127806792cafSBarry Smith + mg - the multigrid context 127906792cafSBarry Smith - n - the number of smoothing steps 128006792cafSBarry Smith 128106792cafSBarry Smith Options Database Key: 1282a2b725a8SWilliam Gropp . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 128306792cafSBarry Smith 128406792cafSBarry Smith Level: advanced 128506792cafSBarry Smith 128695452b02SPatrick Sanan Notes: 128795452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 128806792cafSBarry Smith there is no separate smooth up on the coarsest grid. 128906792cafSBarry Smith 1290710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp() 129106792cafSBarry Smith @*/ 129206792cafSBarry Smith PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 129306792cafSBarry Smith { 129406792cafSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 129506792cafSBarry Smith PC_MG_Levels **mglevels = mg->levels; 129606792cafSBarry Smith PetscErrorCode ierr; 129706792cafSBarry Smith PetscInt i,levels; 129806792cafSBarry Smith 129906792cafSBarry Smith PetscFunctionBegin; 130006792cafSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 130106792cafSBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 13021a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 130306792cafSBarry Smith levels = mglevels[0]->levels; 130406792cafSBarry Smith 130506792cafSBarry Smith for (i=1; i<levels; i++) { 130606792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 130706792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 130806792cafSBarry Smith mg->default_smoothu = n; 130906792cafSBarry Smith mg->default_smoothd = n; 131006792cafSBarry Smith } 131106792cafSBarry Smith PetscFunctionReturn(0); 131206792cafSBarry Smith } 131306792cafSBarry Smith 1314f442ab6aSBarry Smith /*@ 13158e5aa403SBarry Smith PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels 1316710315b6SLawrence Mitchell and adds the suffix _up to the options name 1317f442ab6aSBarry Smith 1318f442ab6aSBarry Smith Logically Collective on PC 1319f442ab6aSBarry Smith 1320f442ab6aSBarry Smith Input Parameters: 1321f442ab6aSBarry Smith . pc - the preconditioner context 1322f442ab6aSBarry Smith 1323f442ab6aSBarry Smith Options Database Key: 1324f442ab6aSBarry Smith . -pc_mg_distinct_smoothup 1325f442ab6aSBarry Smith 1326f442ab6aSBarry Smith Level: advanced 1327f442ab6aSBarry Smith 132895452b02SPatrick Sanan Notes: 132995452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 1330f442ab6aSBarry Smith there is no separate smooth up on the coarsest grid. 1331f442ab6aSBarry Smith 1332710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth() 1333f442ab6aSBarry Smith @*/ 1334f442ab6aSBarry Smith PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1335f442ab6aSBarry Smith { 1336f442ab6aSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1337f442ab6aSBarry Smith PC_MG_Levels **mglevels = mg->levels; 1338f442ab6aSBarry Smith PetscErrorCode ierr; 1339f442ab6aSBarry Smith PetscInt i,levels; 1340f442ab6aSBarry Smith KSP subksp; 1341f442ab6aSBarry Smith 1342f442ab6aSBarry Smith PetscFunctionBegin; 1343f442ab6aSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1344f442ab6aSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1345f442ab6aSBarry Smith levels = mglevels[0]->levels; 1346f442ab6aSBarry Smith 1347f442ab6aSBarry Smith for (i=1; i<levels; i++) { 1348710315b6SLawrence Mitchell const char *prefix = NULL; 1349f442ab6aSBarry Smith /* make sure smoother up and down are different */ 1350f442ab6aSBarry Smith ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr); 1351710315b6SLawrence Mitchell ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr); 1352710315b6SLawrence Mitchell ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr); 1353f442ab6aSBarry Smith ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr); 1354f442ab6aSBarry Smith } 1355f442ab6aSBarry Smith PetscFunctionReturn(0); 1356f442ab6aSBarry Smith } 1357f442ab6aSBarry Smith 135807a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 135907a4832bSFande Kong PetscErrorCode PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[]) 136007a4832bSFande Kong { 136107a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 136207a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 136307a4832bSFande Kong Mat *mat; 136407a4832bSFande Kong PetscInt l; 136507a4832bSFande Kong PetscErrorCode ierr; 136607a4832bSFande Kong 136707a4832bSFande Kong PetscFunctionBegin; 136807a4832bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 136907a4832bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 137007a4832bSFande Kong for (l=1; l< mg->nlevels; l++) { 137107a4832bSFande Kong mat[l-1] = mglevels[l]->interpolate; 137207a4832bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr); 137307a4832bSFande Kong } 137407a4832bSFande Kong *num_levels = mg->nlevels; 137507a4832bSFande Kong *interpolations = mat; 137607a4832bSFande Kong PetscFunctionReturn(0); 137707a4832bSFande Kong } 137807a4832bSFande Kong 137907a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 138007a4832bSFande Kong PetscErrorCode PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[]) 138107a4832bSFande Kong { 138207a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 138307a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 138407a4832bSFande Kong PetscInt l; 138507a4832bSFande Kong Mat *mat; 138607a4832bSFande Kong PetscErrorCode ierr; 138707a4832bSFande Kong 138807a4832bSFande Kong PetscFunctionBegin; 138907a4832bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 139007a4832bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 139107a4832bSFande Kong for (l=0; l<mg->nlevels-1; l++) { 139207a4832bSFande Kong ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr); 139307a4832bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr); 139407a4832bSFande Kong } 139507a4832bSFande Kong *num_levels = mg->nlevels; 139607a4832bSFande Kong *coarseOperators = mat; 139707a4832bSFande Kong PetscFunctionReturn(0); 139807a4832bSFande Kong } 139907a4832bSFande Kong 1400*f3b08a26SMatthew G. Knepley /*@C 1401*f3b08a26SMatthew G. Knepley PCMGRegisterCoarseSpaceConstructor - Adds a method to the PCMG package for coarse space construction. 1402*f3b08a26SMatthew G. Knepley 1403*f3b08a26SMatthew G. Knepley Not collective 1404*f3b08a26SMatthew G. Knepley 1405*f3b08a26SMatthew G. Knepley Input Parameters: 1406*f3b08a26SMatthew G. Knepley + name - name of the constructor 1407*f3b08a26SMatthew G. Knepley - function - constructor routine 1408*f3b08a26SMatthew G. Knepley 1409*f3b08a26SMatthew G. Knepley Notes: 1410*f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1411*f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp) 1412*f3b08a26SMatthew G. Knepley $ pc - The PC object 1413*f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1414*f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1415*f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1416*f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1417*f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1418*f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1419*f3b08a26SMatthew G. Knepley 1420*f3b08a26SMatthew G. Knepley Level: advanced 1421*f3b08a26SMatthew G. Knepley 1422*f3b08a26SMatthew G. Knepley .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister() 1423*f3b08a26SMatthew G. Knepley @*/ 1424*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)) 1425*f3b08a26SMatthew G. Knepley { 1426*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1427*f3b08a26SMatthew G. Knepley 1428*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1429*f3b08a26SMatthew G. Knepley ierr = PCInitializePackage();CHKERRQ(ierr); 1430*f3b08a26SMatthew G. Knepley ierr = PetscFunctionListAdd(&PCMGCoarseList,name,function);CHKERRQ(ierr); 1431*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1432*f3b08a26SMatthew G. Knepley } 1433*f3b08a26SMatthew G. Knepley 1434*f3b08a26SMatthew G. Knepley /*@C 1435*f3b08a26SMatthew G. Knepley PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method. 1436*f3b08a26SMatthew G. Knepley 1437*f3b08a26SMatthew G. Knepley Not collective 1438*f3b08a26SMatthew G. Knepley 1439*f3b08a26SMatthew G. Knepley Input Parameter: 1440*f3b08a26SMatthew G. Knepley . name - name of the constructor 1441*f3b08a26SMatthew G. Knepley 1442*f3b08a26SMatthew G. Knepley Output Parameter: 1443*f3b08a26SMatthew G. Knepley . function - constructor routine 1444*f3b08a26SMatthew G. Knepley 1445*f3b08a26SMatthew G. Knepley Notes: 1446*f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1447*f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp) 1448*f3b08a26SMatthew G. Knepley $ pc - The PC object 1449*f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1450*f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1451*f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1452*f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1453*f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1454*f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1455*f3b08a26SMatthew G. Knepley 1456*f3b08a26SMatthew G. Knepley Level: advanced 1457*f3b08a26SMatthew G. Knepley 1458*f3b08a26SMatthew G. Knepley .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister() 1459*f3b08a26SMatthew G. Knepley @*/ 1460*f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)) 1461*f3b08a26SMatthew G. Knepley { 1462*f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1463*f3b08a26SMatthew G. Knepley 1464*f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1465*f3b08a26SMatthew G. Knepley ierr = PetscFunctionListFind(PCMGCoarseList,name,function);CHKERRQ(ierr); 1466*f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1467*f3b08a26SMatthew G. Knepley } 1468*f3b08a26SMatthew G. Knepley 14694b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 14704b9ad928SBarry Smith 14713b09bd56SBarry Smith /*MC 1472ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 14733b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 14743b09bd56SBarry Smith 14753b09bd56SBarry Smith Options Database Keys: 14763b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1477391689abSStefano Zampini . -pc_mg_cycle_type <v,w> - provide the cycle desired 14788c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 14793b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 1480710315b6SLawrence Mitchell . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 14812134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 14828cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 14838cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1484e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 14858cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 14868cf6037eSBarry Smith to the binary output file called binaryoutput 14873b09bd56SBarry Smith 148895452b02SPatrick Sanan Notes: 14890b3c753eSRichard Tran Mills If one uses a Krylov method such GMRES or CG as the smoother then one must use KSPFGMRES, KSPGCR, or KSPRICHARDSON as the outer Krylov method 14903b09bd56SBarry Smith 14918cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 14928cf6037eSBarry Smith 149323067569SBarry Smith When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 149423067569SBarry Smith is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 149523067569SBarry Smith (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 149623067569SBarry Smith residual is computed at the end of each cycle. 149723067569SBarry Smith 14983b09bd56SBarry Smith Level: intermediate 14993b09bd56SBarry Smith 15008cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1501710315b6SLawrence Mitchell PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), 1502710315b6SLawrence Mitchell PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 150397177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 15040d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 15053b09bd56SBarry Smith M*/ 15063b09bd56SBarry Smith 15078cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 15084b9ad928SBarry Smith { 1509f3fbd535SBarry Smith PC_MG *mg; 1510f3fbd535SBarry Smith PetscErrorCode ierr; 1511f3fbd535SBarry Smith 15124b9ad928SBarry Smith PetscFunctionBegin; 1513b00a9115SJed Brown ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1514f3fbd535SBarry Smith pc->data = (void*)mg; 1515f3fbd535SBarry Smith mg->nlevels = -1; 151610eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 15172134b1e4SBarry Smith mg->galerkin = PC_MG_GALERKIN_NONE; 1518*f3b08a26SMatthew G. Knepley mg->adaptInterpolation = PETSC_FALSE; 1519*f3b08a26SMatthew G. Knepley mg->Nc = -1; 1520*f3b08a26SMatthew G. Knepley mg->eigenvalue = -1; 1521f3fbd535SBarry Smith 152237a44384SMark Adams pc->useAmat = PETSC_TRUE; 152337a44384SMark Adams 15244b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 15254b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1526a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 15274b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 15284b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 15294b9ad928SBarry Smith pc->ops->view = PCView_MG; 1530fb15c04eSBarry Smith 1531*f3b08a26SMatthew G. Knepley ierr = PetscObjectComposedDataRegister(&mg->eigenvalue);CHKERRQ(ierr); 1532fb15c04eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 153397d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr); 153497d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr); 1535fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);CHKERRQ(ierr); 1536fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);CHKERRQ(ierr); 1537*f3b08a26SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);CHKERRQ(ierr); 1538*f3b08a26SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);CHKERRQ(ierr); 15394b9ad928SBarry Smith PetscFunctionReturn(0); 15404b9ad928SBarry Smith } 1541