1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines the multigrid preconditioner interface. 44b9ad928SBarry Smith */ 5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 6*41b6fd38SMatthew G. Knepley #include <petsc/private/kspimpl.h> 71e25c274SJed Brown #include <petscdm.h> 8391689abSStefano Zampini PETSC_INTERN PetscErrorCode PCPreSolveChangeRHS(PC,PetscBool*); 94b9ad928SBarry Smith 10f3b08a26SMatthew G. Knepley /* 11f3b08a26SMatthew G. Knepley Contains the list of registered coarse space construction routines 12f3b08a26SMatthew G. Knepley */ 13f3b08a26SMatthew G. Knepley PetscFunctionList PCMGCoarseList = NULL; 14f3b08a26SMatthew G. Knepley 1531567311SBarry Smith PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) 164b9ad928SBarry Smith { 1731567311SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1831567311SBarry Smith PC_MG_Levels *mgc,*mglevels = *mglevelsin; 196849ba73SBarry Smith PetscErrorCode ierr; 2031567311SBarry Smith PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; 214b9ad928SBarry Smith 224b9ad928SBarry Smith PetscFunctionBegin; 2363e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2431567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 25c0decd05SBarry Smith ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr); 2663e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 2731567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 2863e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 2931567311SBarry Smith ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); 3063e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 314b9ad928SBarry Smith 324b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 3331567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 344b9ad928SBarry Smith PetscReal rnorm; 3531567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 364b9ad928SBarry Smith if (rnorm <= mg->ttol) { 3770441072SBarry Smith if (rnorm < mg->abstol) { 384d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 3957622a8eSBarry 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); 404b9ad928SBarry Smith } else { 414d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 4257622a8eSBarry 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); 434b9ad928SBarry Smith } 444b9ad928SBarry Smith PetscFunctionReturn(0); 454b9ad928SBarry Smith } 464b9ad928SBarry Smith } 474b9ad928SBarry Smith 4831567311SBarry Smith mgc = *(mglevelsin - 1); 4963e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5031567311SBarry Smith ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); 5163e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 52efb30889SBarry Smith ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); 534b9ad928SBarry Smith while (cycles--) { 5431567311SBarry Smith ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); 554b9ad928SBarry Smith } 5663e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5731567311SBarry Smith ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 5863e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5963e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 6031567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 61c0decd05SBarry Smith ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr); 62*41b6fd38SMatthew G. Knepley if (mglevels->cr) { 63*41b6fd38SMatthew G. Knepley /* TODO Turn on copy and turn off noisy if we have an exact solution 64*41b6fd38SMatthew G. Knepley ierr = VecCopy(mglevels->x, mglevels->crx);CHKERRQ(ierr); 65*41b6fd38SMatthew G. Knepley ierr = VecCopy(mglevels->b, mglevels->crb);CHKERRQ(ierr); */ 66*41b6fd38SMatthew G. Knepley ierr = KSPSetNoisy_Private(mglevels->crx);CHKERRQ(ierr); 67*41b6fd38SMatthew G. Knepley ierr = KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx);CHKERRQ(ierr); /* compatible relaxation */ 68*41b6fd38SMatthew G. Knepley ierr = KSPCheckSolve(mglevels->cr,pc,mglevels->crx);CHKERRQ(ierr); 69*41b6fd38SMatthew G. Knepley } 7063e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 714b9ad928SBarry Smith } 724b9ad928SBarry Smith PetscFunctionReturn(0); 734b9ad928SBarry Smith } 744b9ad928SBarry Smith 75ace3abfcSBarry 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) 764b9ad928SBarry Smith { 77f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 78f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 79dfbe8321SBarry Smith PetscErrorCode ierr; 80391689abSStefano Zampini PC tpc; 81391689abSStefano Zampini PetscBool changeu,changed; 82f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 834b9ad928SBarry Smith 844b9ad928SBarry Smith PetscFunctionBegin; 85a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 86a762d673SBarry Smith for (i=0; i<levels; i++) { 87a762d673SBarry Smith if (!mglevels[i]->A) { 88a762d673SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 89a762d673SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 90a762d673SBarry Smith } 91a762d673SBarry Smith } 92391689abSStefano Zampini 93391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 94391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 95391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 96391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 97391689abSStefano Zampini if (!changed && !changeu) { 98391689abSStefano Zampini ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 99f3fbd535SBarry Smith mglevels[levels-1]->b = b; 100391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 101391689abSStefano Zampini if (!mglevels[levels-1]->b) { 102391689abSStefano Zampini Vec *vec; 103391689abSStefano Zampini 104391689abSStefano Zampini ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 105391689abSStefano Zampini mglevels[levels-1]->b = *vec; 1067ae21d43SStefano Zampini ierr = PetscFree(vec);CHKERRQ(ierr); 107391689abSStefano Zampini } 108391689abSStefano Zampini ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 109391689abSStefano Zampini } 110f3fbd535SBarry Smith mglevels[levels-1]->x = x; 1114b9ad928SBarry Smith 11231567311SBarry Smith mg->rtol = rtol; 11331567311SBarry Smith mg->abstol = abstol; 11431567311SBarry Smith mg->dtol = dtol; 1154b9ad928SBarry Smith if (rtol) { 1164b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 1174b9ad928SBarry Smith PetscReal rnorm; 1187319c654SBarry Smith if (zeroguess) { 1197319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 1207319c654SBarry Smith } else { 121f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 1224b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 1237319c654SBarry Smith } 12431567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 1252fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 1262fa5cd67SKarl Rupp else mg->ttol = 0.0; 1274b9ad928SBarry Smith 1284d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 129336babb1SJed Brown stop prematurely due to small residual */ 1304d0a8057SBarry Smith for (i=1; i<levels; i++) { 131f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 132f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 13323067569SBarry Smith /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 13423067569SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 135f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 1364b9ad928SBarry Smith } 1374d0a8057SBarry Smith } 1384d0a8057SBarry Smith 1394d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1404d0a8057SBarry Smith for (i=0; i<its; i++) { 141f3fbd535SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); 1424d0a8057SBarry Smith if (*reason) break; 1434d0a8057SBarry Smith } 1444d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1454d0a8057SBarry Smith *outits = i; 146391689abSStefano Zampini if (!changed && !changeu) mglevels[levels-1]->b = NULL; 1474b9ad928SBarry Smith PetscFunctionReturn(0); 1484b9ad928SBarry Smith } 1494b9ad928SBarry Smith 1503aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 1513aeaf226SBarry Smith { 1523aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1533aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 1543aeaf226SBarry Smith PetscErrorCode ierr; 155f3b08a26SMatthew G. Knepley PetscInt i,c,n; 1563aeaf226SBarry Smith 1573aeaf226SBarry Smith PetscFunctionBegin; 1583aeaf226SBarry Smith if (mglevels) { 1593aeaf226SBarry Smith n = mglevels[0]->levels; 1603aeaf226SBarry Smith for (i=0; i<n-1; i++) { 1613aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 1623aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 1633aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 164*41b6fd38SMatthew G. Knepley ierr = VecDestroy(&mglevels[i]->crx);CHKERRQ(ierr); 165*41b6fd38SMatthew G. Knepley ierr = VecDestroy(&mglevels[i]->crb);CHKERRQ(ierr); 1663aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 1673aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 1680de8493bSLawrence Mitchell ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr); 16973250ac0SBarry Smith ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 1703aeaf226SBarry Smith } 171*41b6fd38SMatthew G. Knepley ierr = VecDestroy(&mglevels[n-1]->crx);CHKERRQ(ierr); 172*41b6fd38SMatthew G. Knepley ierr = VecDestroy(&mglevels[n-1]->crb);CHKERRQ(ierr); 173391689abSStefano Zampini /* this is not null only if the smoother on the finest level 174391689abSStefano Zampini changes the rhs during PreSolve */ 175391689abSStefano Zampini ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr); 1763aeaf226SBarry Smith 1773aeaf226SBarry Smith for (i=0; i<n; i++) { 178f3b08a26SMatthew G. Knepley if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {ierr = VecDestroy(&mglevels[i]->coarseSpace[c]);CHKERRQ(ierr);} 179f3b08a26SMatthew G. Knepley ierr = PetscFree(mglevels[i]->coarseSpace);CHKERRQ(ierr); 180f3b08a26SMatthew G. Knepley mglevels[i]->coarseSpace = NULL; 1813aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 1823aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1833aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 1843aeaf226SBarry Smith } 1853aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 186*41b6fd38SMatthew G. Knepley if (mglevels[i]->cr) {ierr = KSPReset(mglevels[i]->cr);CHKERRQ(ierr);} 1873aeaf226SBarry Smith } 188f3b08a26SMatthew G. Knepley mg->Nc = 0; 1893aeaf226SBarry Smith } 1903aeaf226SBarry Smith PetscFunctionReturn(0); 1913aeaf226SBarry Smith } 1923aeaf226SBarry Smith 193*41b6fd38SMatthew G. Knepley /* Implementing CR 194*41b6fd38SMatthew G. Knepley 195*41b6fd38SMatthew G. Knepley We only want to make corrections that ``do not change'' the coarse solution. What we mean by not changing is that if I prolong my coarse solution to the fine grid and then inject that fine solution back to the coarse grid, I get the same answer. Injection is what Brannick calls R. We want the complementary projector to Inj, which we will call S, after Brannick, so that Inj S = 0. Now the orthogonal projector onto the range of Inj^T is 196*41b6fd38SMatthew G. Knepley 197*41b6fd38SMatthew G. Knepley Inj^T (Inj Inj^T)^{-1} Inj 198*41b6fd38SMatthew G. Knepley 199*41b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have 200*41b6fd38SMatthew G. Knepley 201*41b6fd38SMatthew G. Knepley Inj^T Inj 202*41b6fd38SMatthew G. Knepley 203*41b6fd38SMatthew G. Knepley and 204*41b6fd38SMatthew G. Knepley 205*41b6fd38SMatthew G. Knepley S = I - Inj^T Inj 206*41b6fd38SMatthew G. Knepley 207*41b6fd38SMatthew G. Knepley since 208*41b6fd38SMatthew G. Knepley 209*41b6fd38SMatthew G. Knepley Inj S = Inj - (Inj Inj^T) Inj = 0. 210*41b6fd38SMatthew G. Knepley 211*41b6fd38SMatthew G. Knepley Brannick suggests 212*41b6fd38SMatthew G. Knepley 213*41b6fd38SMatthew G. Knepley A \to S^T A S \qquad\mathrm{and}\qquad M \to S^T M S 214*41b6fd38SMatthew G. Knepley 215*41b6fd38SMatthew G. Knepley but I do not think his :math:`S^T S = I` is correct. Our S is an orthogonal projector, so :math:`S^T S = S^2 = S`. We will use 216*41b6fd38SMatthew G. Knepley 217*41b6fd38SMatthew G. Knepley M^{-1} A \to S M^{-1} A S 218*41b6fd38SMatthew G. Knepley 219*41b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left. 220*41b6fd38SMatthew G. Knepley 221*41b6fd38SMatthew G. Knepley Check: || Inj P - I ||_F < tol 222*41b6fd38SMatthew G. Knepley Check: In general, Inj Inj^T = I 223*41b6fd38SMatthew G. Knepley */ 224*41b6fd38SMatthew G. Knepley 225*41b6fd38SMatthew G. Knepley typedef struct { 226*41b6fd38SMatthew G. Knepley PC mg; /* The PCMG object */ 227*41b6fd38SMatthew G. Knepley PetscInt l; /* The multigrid level for this solver */ 228*41b6fd38SMatthew G. Knepley Mat Inj; /* The injection matrix */ 229*41b6fd38SMatthew G. Knepley Mat S; /* I - Inj^T Inj */ 230*41b6fd38SMatthew G. Knepley } CRContext; 231*41b6fd38SMatthew G. Knepley 232*41b6fd38SMatthew G. Knepley static PetscErrorCode CRSetup_Private(PC pc) 233*41b6fd38SMatthew G. Knepley { 234*41b6fd38SMatthew G. Knepley CRContext *ctx; 235*41b6fd38SMatthew G. Knepley Mat It; 236*41b6fd38SMatthew G. Knepley PetscErrorCode ierr; 237*41b6fd38SMatthew G. Knepley 238*41b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 239*41b6fd38SMatthew G. Knepley ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr); 240*41b6fd38SMatthew G. Knepley ierr = PCMGGetInjection(ctx->mg, ctx->l, &It);CHKERRQ(ierr); 241*41b6fd38SMatthew G. Knepley if (!It) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG"); 242*41b6fd38SMatthew G. Knepley ierr = MatCreateTranspose(It, &ctx->Inj);CHKERRQ(ierr); 243*41b6fd38SMatthew G. Knepley ierr = MatCreateNormal(ctx->Inj, &ctx->S);CHKERRQ(ierr); 244*41b6fd38SMatthew G. Knepley ierr = MatScale(ctx->S, -1.0);CHKERRQ(ierr); 245*41b6fd38SMatthew G. Knepley ierr = MatShift(ctx->S, 1.0);CHKERRQ(ierr); 246*41b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 247*41b6fd38SMatthew G. Knepley } 248*41b6fd38SMatthew G. Knepley 249*41b6fd38SMatthew G. Knepley static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y) 250*41b6fd38SMatthew G. Knepley { 251*41b6fd38SMatthew G. Knepley CRContext *ctx; 252*41b6fd38SMatthew G. Knepley PetscErrorCode ierr; 253*41b6fd38SMatthew G. Knepley 254*41b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 255*41b6fd38SMatthew G. Knepley ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr); 256*41b6fd38SMatthew G. Knepley ierr = MatMult(ctx->S, x, y);CHKERRQ(ierr); 257*41b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 258*41b6fd38SMatthew G. Knepley } 259*41b6fd38SMatthew G. Knepley 260*41b6fd38SMatthew G. Knepley static PetscErrorCode CRDestroy_Private(PC pc) 261*41b6fd38SMatthew G. Knepley { 262*41b6fd38SMatthew G. Knepley CRContext *ctx; 263*41b6fd38SMatthew G. Knepley PetscErrorCode ierr; 264*41b6fd38SMatthew G. Knepley 265*41b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 266*41b6fd38SMatthew G. Knepley ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr); 267*41b6fd38SMatthew G. Knepley ierr = MatDestroy(&ctx->Inj);CHKERRQ(ierr); 268*41b6fd38SMatthew G. Knepley ierr = MatDestroy(&ctx->S);CHKERRQ(ierr); 269*41b6fd38SMatthew G. Knepley ierr = PetscFree(ctx);CHKERRQ(ierr); 270*41b6fd38SMatthew G. Knepley ierr = PCShellSetContext(pc, NULL);CHKERRQ(ierr); 271*41b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 272*41b6fd38SMatthew G. Knepley } 273*41b6fd38SMatthew G. Knepley 274*41b6fd38SMatthew G. Knepley static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr) 275*41b6fd38SMatthew G. Knepley { 276*41b6fd38SMatthew G. Knepley CRContext *ctx; 277*41b6fd38SMatthew G. Knepley PetscErrorCode ierr; 278*41b6fd38SMatthew G. Knepley 279*41b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 280*41b6fd38SMatthew G. Knepley ierr = PCCreate(PetscObjectComm((PetscObject) pc), cr);CHKERRQ(ierr); 281*41b6fd38SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)");CHKERRQ(ierr); 282*41b6fd38SMatthew G. Knepley ierr = PetscCalloc1(1, &ctx);CHKERRQ(ierr); 283*41b6fd38SMatthew G. Knepley ctx->mg = pc; 284*41b6fd38SMatthew G. Knepley ctx->l = l; 285*41b6fd38SMatthew G. Knepley ierr = PCSetType(*cr, PCSHELL);CHKERRQ(ierr); 286*41b6fd38SMatthew G. Knepley ierr = PCShellSetContext(*cr, ctx);CHKERRQ(ierr); 287*41b6fd38SMatthew G. Knepley ierr = PCShellSetApply(*cr, CRApply_Private);CHKERRQ(ierr); 288*41b6fd38SMatthew G. Knepley ierr = PCShellSetSetUp(*cr, CRSetup_Private);CHKERRQ(ierr); 289*41b6fd38SMatthew G. Knepley ierr = PCShellSetDestroy(*cr, CRDestroy_Private);CHKERRQ(ierr); 290*41b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 291*41b6fd38SMatthew G. Knepley } 292*41b6fd38SMatthew G. Knepley 29397d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms) 2944b9ad928SBarry Smith { 295dfbe8321SBarry Smith PetscErrorCode ierr; 296f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 297ce94432eSBarry Smith MPI_Comm comm; 2983aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 29910eca3edSLisandro Dalcin PCMGType mgtype = mg->am; 30010167fecSBarry Smith PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V; 301f3fbd535SBarry Smith PetscInt i; 302f3fbd535SBarry Smith PetscMPIInt size; 303f3fbd535SBarry Smith const char *prefix; 304f3fbd535SBarry Smith PC ipc; 3053aeaf226SBarry Smith PetscInt n; 3064b9ad928SBarry Smith 3074b9ad928SBarry Smith PetscFunctionBegin; 3080700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 309548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 310548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 3111a97d7b6SStefano Zampini ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 3123aeaf226SBarry Smith if (mglevels) { 31310eca3edSLisandro Dalcin mgctype = mglevels[0]->cycles; 3143aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 3153aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 3163aeaf226SBarry Smith n = mglevels[0]->levels; 3173aeaf226SBarry Smith for (i=0; i<n; i++) { 3183aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 3193aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 3203aeaf226SBarry Smith } 3213aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 322*41b6fd38SMatthew G. Knepley ierr = KSPDestroy(&mglevels[i]->cr);CHKERRQ(ierr); 3233aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 3243aeaf226SBarry Smith } 3253aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 3263aeaf226SBarry Smith } 327f3fbd535SBarry Smith 328f3fbd535SBarry Smith mg->nlevels = levels; 329f3fbd535SBarry Smith 330785e854fSJed Brown ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr); 3313bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 332f3fbd535SBarry Smith 333f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 334f3fbd535SBarry Smith 335a9db87a2SMatthew G Knepley mg->stageApply = 0; 336f3fbd535SBarry Smith for (i=0; i<levels; i++) { 337b00a9115SJed Brown ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr); 3382fa5cd67SKarl Rupp 339f3fbd535SBarry Smith mglevels[i]->level = i; 340f3fbd535SBarry Smith mglevels[i]->levels = levels; 34110eca3edSLisandro Dalcin mglevels[i]->cycles = mgctype; 342336babb1SJed Brown mg->default_smoothu = 2; 343336babb1SJed Brown mg->default_smoothd = 2; 34463e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 34563e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 34663e6d426SJed Brown mglevels[i]->eventresidual = 0; 34763e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 348f3fbd535SBarry Smith 349f3fbd535SBarry Smith if (comms) comm = comms[i]; 350f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 351422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr); 3520ee9a668SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 3530ee9a668SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 3540ee9a668SBarry Smith ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 3550ee9a668SBarry Smith if (i || levels == 1) { 3560ee9a668SBarry Smith char tprefix[128]; 3570ee9a668SBarry Smith 358336babb1SJed Brown ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 3590059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 360336babb1SJed Brown ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 361336babb1SJed Brown ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 362336babb1SJed Brown ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 3630ee9a668SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 364f3fbd535SBarry Smith 365*41b6fd38SMatthew G. Knepley ierr = PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);CHKERRQ(ierr); 3660ee9a668SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 3670ee9a668SBarry Smith } else { 368f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 369f3fbd535SBarry Smith 3709dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 371f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 372f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 373f3fbd535SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 374f3fbd535SBarry Smith if (size > 1) { 375f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 376f3fbd535SBarry Smith } else { 377f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 378f3fbd535SBarry Smith } 379753b7fb9SBarry Smith ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 380f3fbd535SBarry Smith } 3813bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr); 3822fa5cd67SKarl Rupp 383f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 38431567311SBarry Smith mg->rtol = 0.0; 38531567311SBarry Smith mg->abstol = 0.0; 38631567311SBarry Smith mg->dtol = 0.0; 38731567311SBarry Smith mg->ttol = 0.0; 38831567311SBarry Smith mg->cyclesperpcapply = 1; 389f3fbd535SBarry Smith } 390f3fbd535SBarry Smith mg->levels = mglevels; 39110eca3edSLisandro Dalcin ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 3924b9ad928SBarry Smith PetscFunctionReturn(0); 3934b9ad928SBarry Smith } 3944b9ad928SBarry Smith 39597d33e41SMatthew G. Knepley /*@C 39697d33e41SMatthew G. Knepley PCMGSetLevels - Sets the number of levels to use with MG. 39797d33e41SMatthew G. Knepley Must be called before any other MG routine. 39897d33e41SMatthew G. Knepley 39997d33e41SMatthew G. Knepley Logically Collective on PC 40097d33e41SMatthew G. Knepley 40197d33e41SMatthew G. Knepley Input Parameters: 40297d33e41SMatthew G. Knepley + pc - the preconditioner context 40397d33e41SMatthew G. Knepley . levels - the number of levels 40497d33e41SMatthew G. Knepley - comms - optional communicators for each level; this is to allow solving the coarser problems 40597d33e41SMatthew G. Knepley on smaller sets of processors. 40697d33e41SMatthew G. Knepley 40797d33e41SMatthew G. Knepley Level: intermediate 40897d33e41SMatthew G. Knepley 40997d33e41SMatthew G. Knepley Notes: 41097d33e41SMatthew G. Knepley If the number of levels is one then the multigrid uses the -mg_levels prefix 41197d33e41SMatthew G. Knepley for setting the level options rather than the -mg_coarse prefix. 41297d33e41SMatthew G. Knepley 41397d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels() 41497d33e41SMatthew G. Knepley @*/ 41597d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 41697d33e41SMatthew G. Knepley { 41797d33e41SMatthew G. Knepley PetscErrorCode ierr; 41897d33e41SMatthew G. Knepley 41997d33e41SMatthew G. Knepley PetscFunctionBegin; 42097d33e41SMatthew G. Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 42197d33e41SMatthew G. Knepley if (comms) PetscValidPointer(comms,3); 42237b5128cSMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr); 42397d33e41SMatthew G. Knepley PetscFunctionReturn(0); 42497d33e41SMatthew G. Knepley } 42597d33e41SMatthew G. Knepley 426c07bf074SBarry Smith 427c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 428c07bf074SBarry Smith { 429c07bf074SBarry Smith PetscErrorCode ierr; 430a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 431a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 432a06653b4SBarry Smith PetscInt i,n; 433c07bf074SBarry Smith 434c07bf074SBarry Smith PetscFunctionBegin; 435a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 436a06653b4SBarry Smith if (mglevels) { 437a06653b4SBarry Smith n = mglevels[0]->levels; 438a06653b4SBarry Smith for (i=0; i<n; i++) { 439a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 4406bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 441a06653b4SBarry Smith } 4426bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 443*41b6fd38SMatthew G. Knepley ierr = KSPDestroy(&mglevels[i]->cr);CHKERRQ(ierr); 444a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 445a06653b4SBarry Smith } 446a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 447a06653b4SBarry Smith } 448c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 449fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);CHKERRQ(ierr); 450fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);CHKERRQ(ierr); 451f3fbd535SBarry Smith PetscFunctionReturn(0); 452f3fbd535SBarry Smith } 453f3fbd535SBarry Smith 454f3fbd535SBarry Smith 455f3fbd535SBarry Smith 45609573ac7SBarry Smith extern PetscErrorCode PCMGACycle_Private(PC,PC_MG_Levels**); 45709573ac7SBarry Smith extern PetscErrorCode PCMGFCycle_Private(PC,PC_MG_Levels**); 45809573ac7SBarry Smith extern PetscErrorCode PCMGKCycle_Private(PC,PC_MG_Levels**); 459f3fbd535SBarry Smith 460f3fbd535SBarry Smith /* 461f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 462f3fbd535SBarry Smith or full cycle of multigrid. 463f3fbd535SBarry Smith 464f3fbd535SBarry Smith Note: 465f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 466f3fbd535SBarry Smith */ 467f3fbd535SBarry Smith static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 468f3fbd535SBarry Smith { 469f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 470f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 471f3fbd535SBarry Smith PetscErrorCode ierr; 472391689abSStefano Zampini PC tpc; 473f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 474391689abSStefano Zampini PetscBool changeu,changed; 475f3fbd535SBarry Smith 476f3fbd535SBarry Smith PetscFunctionBegin; 477a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 478e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 479298cc208SBarry Smith for (i=0; i<levels; i++) { 480e1d8e5deSBarry Smith if (!mglevels[i]->A) { 48123ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 482298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 483e1d8e5deSBarry Smith } 484e1d8e5deSBarry Smith } 485e1d8e5deSBarry Smith 486391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 487391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 488391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 489391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 490391689abSStefano Zampini if (!changeu && !changed) { 491391689abSStefano Zampini ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 492f3fbd535SBarry Smith mglevels[levels-1]->b = b; 493391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 494391689abSStefano Zampini if (!mglevels[levels-1]->b) { 495391689abSStefano Zampini Vec *vec; 496391689abSStefano Zampini 497391689abSStefano Zampini ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 498391689abSStefano Zampini mglevels[levels-1]->b = *vec; 4997ae21d43SStefano Zampini ierr = PetscFree(vec);CHKERRQ(ierr); 500391689abSStefano Zampini } 501391689abSStefano Zampini ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 502391689abSStefano Zampini } 503f3fbd535SBarry Smith mglevels[levels-1]->x = x; 504391689abSStefano Zampini 50531567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 506f3fbd535SBarry Smith ierr = VecSet(x,0.0);CHKERRQ(ierr); 50731567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 5080298fd71SBarry Smith ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr); 509f3fbd535SBarry Smith } 5102fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 51131567311SBarry Smith ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); 5122fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 51331567311SBarry Smith ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); 5142fa5cd67SKarl Rupp } else { 515f3fbd535SBarry Smith ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); 516f3fbd535SBarry Smith } 517a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 518391689abSStefano Zampini if (!changeu && !changed) mglevels[levels-1]->b = NULL; 519f3fbd535SBarry Smith PetscFunctionReturn(0); 520f3fbd535SBarry Smith } 521f3fbd535SBarry Smith 522f3fbd535SBarry Smith 5234416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 524f3fbd535SBarry Smith { 525f3fbd535SBarry Smith PetscErrorCode ierr; 526710315b6SLawrence Mitchell PetscInt levels,cycles; 527f3b08a26SMatthew G. Knepley PetscBool flg, flg2; 528f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 5293d3eaba7SBarry Smith PC_MG_Levels **mglevels; 530f3fbd535SBarry Smith PCMGType mgtype; 531f3fbd535SBarry Smith PCMGCycleType mgctype; 5322134b1e4SBarry Smith PCMGGalerkinType gtype; 533f3fbd535SBarry Smith 534f3fbd535SBarry Smith PetscFunctionBegin; 5351d743356SStefano Zampini levels = PetscMax(mg->nlevels,1); 536e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr); 537f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 5381a97d7b6SStefano Zampini if (!flg && !mg->levels && pc->dm) { 539eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 540eb3f98d2SBarry Smith levels++; 5413aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 542eb3f98d2SBarry Smith } 5430298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 5443aeaf226SBarry Smith mglevels = mg->levels; 5453aeaf226SBarry Smith 546f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 547f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 548f3fbd535SBarry Smith if (flg) { 549f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 5502fa5cd67SKarl Rupp } 5512134b1e4SBarry Smith gtype = mg->galerkin; 5522134b1e4SBarry Smith ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg);CHKERRQ(ierr); 5532134b1e4SBarry Smith if (flg) { 5542134b1e4SBarry Smith ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr); 555f3fbd535SBarry Smith } 556f3b08a26SMatthew G. Knepley flg2 = PETSC_FALSE; 557f3b08a26SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr); 558f3b08a26SMatthew G. Knepley if (flg) {ierr = PCMGSetAdaptInterpolation(pc, flg2);CHKERRQ(ierr);} 559f3b08a26SMatthew 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); 560f3b08a26SMatthew 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); 561f3b08a26SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);CHKERRQ(ierr); 562*41b6fd38SMatthew G. Knepley flg2 = PETSC_FALSE; 563*41b6fd38SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr); 564*41b6fd38SMatthew G. Knepley if (flg) {ierr = PCMGSetAdaptCR(pc, flg2);CHKERRQ(ierr);} 565f56b1265SBarry Smith flg = PETSC_FALSE; 5668e5aa403SBarry Smith ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr); 567f442ab6aSBarry Smith if (flg) { 568f442ab6aSBarry Smith ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr); 569f442ab6aSBarry Smith } 57031567311SBarry Smith mgtype = mg->am; 571f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 572f3fbd535SBarry Smith if (flg) { 573f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 574f3fbd535SBarry Smith } 57531567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 5765363c1e0SLisandro Dalcin ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 577f3fbd535SBarry Smith if (flg) { 578f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 579f3fbd535SBarry Smith } 580f3fbd535SBarry Smith } 581f3fbd535SBarry Smith flg = PETSC_FALSE; 5820298fd71SBarry Smith ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 583f3fbd535SBarry Smith if (flg) { 584f3fbd535SBarry Smith PetscInt i; 585f3fbd535SBarry Smith char eventname[128]; 5861a97d7b6SStefano Zampini 587f3fbd535SBarry Smith levels = mglevels[0]->levels; 588f3fbd535SBarry Smith for (i=0; i<levels; i++) { 589f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 59063e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 591f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 59263e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 593f3fbd535SBarry Smith if (i) { 594f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 59563e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 596f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 59763e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 598f3fbd535SBarry Smith } 599f3fbd535SBarry Smith } 600eec35531SMatthew G Knepley 601ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 602eec35531SMatthew G Knepley { 603eec35531SMatthew G Knepley const char *sname = "MG Apply"; 604eec35531SMatthew G Knepley PetscStageLog stageLog; 605eec35531SMatthew G Knepley PetscInt st; 606eec35531SMatthew G Knepley 607eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 608eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 609eec35531SMatthew G Knepley PetscBool same; 610eec35531SMatthew G Knepley 611eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 6122fa5cd67SKarl Rupp if (same) mg->stageApply = st; 613eec35531SMatthew G Knepley } 614eec35531SMatthew G Knepley if (!mg->stageApply) { 615eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 616eec35531SMatthew G Knepley } 617eec35531SMatthew G Knepley } 618ec60ca73SSatish Balay #endif 619f3fbd535SBarry Smith } 620f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 621f3b08a26SMatthew G. Knepley /* Check option consistency */ 622f3b08a26SMatthew G. Knepley ierr = PCMGGetGalerkin(pc, >ype);CHKERRQ(ierr); 623f3b08a26SMatthew G. Knepley ierr = PCMGGetAdaptInterpolation(pc, &flg);CHKERRQ(ierr); 624f3b08a26SMatthew 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"); 625f3fbd535SBarry Smith PetscFunctionReturn(0); 626f3fbd535SBarry Smith } 627f3fbd535SBarry Smith 6280a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL}; 6290a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL}; 6300a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL}; 631f3b08a26SMatthew G. Knepley const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL}; 632f3fbd535SBarry Smith 6339804daf3SBarry Smith #include <petscdraw.h> 634f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 635f3fbd535SBarry Smith { 636f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 637f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 638f3fbd535SBarry Smith PetscErrorCode ierr; 639e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 640219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 641f3fbd535SBarry Smith 642f3fbd535SBarry Smith PetscFunctionBegin; 643251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 6445b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 645219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 646f3fbd535SBarry Smith if (iascii) { 647e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 648efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 64931567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 65031567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 651f3fbd535SBarry Smith } 6522134b1e4SBarry Smith if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 653f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 6542134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 6552134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr); 6562134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 6572134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr); 6582134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 6592134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr); 6604f66f45eSBarry Smith } else { 6614f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 662f3fbd535SBarry Smith } 6635adeb434SBarry Smith if (mg->view){ 6645adeb434SBarry Smith ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr); 6655adeb434SBarry Smith } 666f3fbd535SBarry Smith for (i=0; i<levels; i++) { 667f3fbd535SBarry Smith if (!i) { 66889cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 669f3fbd535SBarry Smith } else { 67089cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 671f3fbd535SBarry Smith } 672f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 673f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 674f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 675f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 676f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 677f3fbd535SBarry Smith } else if (i) { 67889cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 679f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 680f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 681f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 682f3fbd535SBarry Smith } 683*41b6fd38SMatthew G. Knepley if (i && mglevels[i]->cr) { 684*41b6fd38SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer,"CR solver on level %D -------------------------------\n",i);CHKERRQ(ierr); 685*41b6fd38SMatthew G. Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 686*41b6fd38SMatthew G. Knepley ierr = KSPView(mglevels[i]->cr,viewer);CHKERRQ(ierr); 687*41b6fd38SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 688*41b6fd38SMatthew G. Knepley } 689f3fbd535SBarry Smith } 6905b0b0462SBarry Smith } else if (isbinary) { 6915b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 6925b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 6935b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 6945b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 6955b0b0462SBarry Smith } 6965b0b0462SBarry Smith } 697219b1045SBarry Smith } else if (isdraw) { 698219b1045SBarry Smith PetscDraw draw; 699b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 700219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 701219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 7020298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 703219b1045SBarry Smith bottom = y - th; 704219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 705b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 706219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 707219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 708219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 709b4375e8dSPeter Brune } else { 710b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 711b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 712b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 713b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 714b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 715b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 716b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 717b4375e8dSPeter Brune } 7180298fd71SBarry Smith ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 7191cd381d6SBarry Smith bottom -= th; 720219b1045SBarry Smith } 7215b0b0462SBarry Smith } 722f3fbd535SBarry Smith PetscFunctionReturn(0); 723f3fbd535SBarry Smith } 724f3fbd535SBarry Smith 725af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 726cab2e9ccSBarry Smith 727f3fbd535SBarry Smith /* 728f3fbd535SBarry Smith Calls setup for the KSP on each level 729f3fbd535SBarry Smith */ 730f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 731f3fbd535SBarry Smith { 732f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 733f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 734f3fbd535SBarry Smith PetscErrorCode ierr; 7351a97d7b6SStefano Zampini PetscInt i,n; 73698e04984SBarry Smith PC cpc; 7373db492dfSBarry Smith PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 738f3fbd535SBarry Smith Mat dA,dB; 739f3fbd535SBarry Smith Vec tvec; 740218a07d4SBarry Smith DM *dms; 7410a545947SLisandro Dalcin PetscViewer viewer = NULL; 742*41b6fd38SMatthew G. Knepley PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE; 743f3fbd535SBarry Smith 744f3fbd535SBarry Smith PetscFunctionBegin; 7451a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up"); 7461a97d7b6SStefano Zampini n = mglevels[0]->levels; 74701bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 7483aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 7493aeaf226SBarry Smith PetscInt levels; 7503aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 7513aeaf226SBarry Smith levels++; 7523aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 7530298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 7543aeaf226SBarry Smith n = levels; 7553aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 7563aeaf226SBarry Smith mglevels = mg->levels; 7573aeaf226SBarry Smith } 7583aeaf226SBarry Smith } 75998e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 7603aeaf226SBarry Smith 761f3fbd535SBarry Smith 762f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 763f3fbd535SBarry Smith /* so use those from global PC */ 764f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 7650298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 76698e04984SBarry Smith if (opsset) { 76798e04984SBarry Smith Mat mmat; 76823ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 76998e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 77098e04984SBarry Smith } 7714a5f07a5SMark F. Adams 772*41b6fd38SMatthew G. Knepley /* Create CR solvers */ 773*41b6fd38SMatthew G. Knepley ierr = PCMGGetAdaptCR(pc, &doCR);CHKERRQ(ierr); 774*41b6fd38SMatthew G. Knepley if (doCR) { 775*41b6fd38SMatthew G. Knepley const char *prefix; 776*41b6fd38SMatthew G. Knepley 777*41b6fd38SMatthew G. Knepley ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 778*41b6fd38SMatthew G. Knepley for (i = 1; i < n; ++i) { 779*41b6fd38SMatthew G. Knepley PC ipc, cr; 780*41b6fd38SMatthew G. Knepley char crprefix[128]; 781*41b6fd38SMatthew G. Knepley 782*41b6fd38SMatthew G. Knepley ierr = KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr);CHKERRQ(ierr); 783*41b6fd38SMatthew G. Knepley ierr = KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE);CHKERRQ(ierr); 784*41b6fd38SMatthew G. Knepley ierr = PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i);CHKERRQ(ierr); 785*41b6fd38SMatthew G. Knepley ierr = KSPSetOptionsPrefix(mglevels[i]->cr, prefix);CHKERRQ(ierr); 786*41b6fd38SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 787*41b6fd38SMatthew G. Knepley ierr = KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV);CHKERRQ(ierr); 788*41b6fd38SMatthew G. Knepley ierr = KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL);CHKERRQ(ierr); 789*41b6fd38SMatthew G. Knepley ierr = KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED);CHKERRQ(ierr); 790*41b6fd38SMatthew G. Knepley ierr = KSPGetPC(mglevels[i]->cr, &ipc);CHKERRQ(ierr); 791*41b6fd38SMatthew G. Knepley 792*41b6fd38SMatthew G. Knepley ierr = PCSetType(ipc, PCCOMPOSITE);CHKERRQ(ierr); 793*41b6fd38SMatthew G. Knepley ierr = PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 794*41b6fd38SMatthew G. Knepley ierr = PCCompositeAddPCType(ipc, PCSOR);CHKERRQ(ierr); 795*41b6fd38SMatthew G. Knepley ierr = CreateCR_Private(pc, i, &cr);CHKERRQ(ierr); 796*41b6fd38SMatthew G. Knepley ierr = PCCompositeAddPC(ipc, cr);CHKERRQ(ierr); 797*41b6fd38SMatthew G. Knepley ierr = PCDestroy(&cr);CHKERRQ(ierr); 798*41b6fd38SMatthew G. Knepley 799*41b6fd38SMatthew G. Knepley ierr = KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 800*41b6fd38SMatthew G. Knepley ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE);CHKERRQ(ierr); 801*41b6fd38SMatthew G. Knepley ierr = PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i);CHKERRQ(ierr); 802*41b6fd38SMatthew G. Knepley ierr = KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix);CHKERRQ(ierr); 803*41b6fd38SMatthew G. Knepley ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr);CHKERRQ(ierr); 804*41b6fd38SMatthew G. Knepley } 805*41b6fd38SMatthew G. Knepley } 806*41b6fd38SMatthew G. Knepley 8074a5f07a5SMark F. Adams if (!opsset) { 80871b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 8094a5f07a5SMark F. Adams if (use_amat) { 810f3fbd535SBarry 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); 81123ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 81269858f1bSStefano Zampini } else { 8134a5f07a5SMark 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); 81423ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 8154a5f07a5SMark F. Adams } 8164a5f07a5SMark F. Adams } 817f3fbd535SBarry Smith 8189c7ffb39SBarry Smith for (i=n-1; i>0; i--) { 8199c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 8209c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 8219c7ffb39SBarry Smith continue; 8229c7ffb39SBarry Smith } 8239c7ffb39SBarry Smith } 8242134b1e4SBarry Smith 8252134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 8262134b1e4SBarry Smith if (dA == dB) dAeqdB = PETSC_TRUE; 8272134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) { 8282134b1e4SBarry Smith needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 8292134b1e4SBarry Smith } 8302134b1e4SBarry Smith 8312134b1e4SBarry Smith 8329c7ffb39SBarry Smith /* 8339c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 8349c7ffb39SBarry Smith Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 8359c7ffb39SBarry Smith */ 8362134b1e4SBarry Smith if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 8372d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 8382e499ae9SBarry Smith Mat p; 83973250ac0SBarry Smith Vec rscale; 840785e854fSJed Brown ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 841218a07d4SBarry Smith dms[n-1] = pc->dm; 842ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 843ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 84441fce8fdSHong Zhang /* 84541fce8fdSHong Zhang Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level. 84641fce8fdSHong Zhang Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options. 84741fce8fdSHong Zhang But it is safe to use -dm_mat_type. 84841fce8fdSHong Zhang 84941fce8fdSHong Zhang The mat type should not be hardcoded like this, we need to find a better way. 85041fce8fdSHong Zhang ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr); 85141fce8fdSHong Zhang */ 852218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 853942e3340SBarry Smith DMKSP kdm; 854eab52d2bSLawrence Mitchell PetscBool dmhasrestrict, dmhasinject; 8553c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 8562134b1e4SBarry Smith if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 857c27ee7a3SPatrick Farrell if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 858c27ee7a3SPatrick Farrell ierr = KSPSetDM(mglevels[i]->smoothu,dms[i]);CHKERRQ(ierr); 859c27ee7a3SPatrick Farrell if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);CHKERRQ(ierr);} 860c27ee7a3SPatrick Farrell } 861*41b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 862*41b6fd38SMatthew G. Knepley ierr = KSPSetDM(mglevels[i]->cr,dms[i]);CHKERRQ(ierr); 863*41b6fd38SMatthew G. Knepley if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE);CHKERRQ(ierr);} 864*41b6fd38SMatthew G. Knepley } 865942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 866d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 867d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 8680298fd71SBarry Smith kdm->ops->computerhs = NULL; 8690298fd71SBarry Smith kdm->rhsctx = NULL; 87024c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 871e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 8722d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 87300ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 87473250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 8756bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 876218a07d4SBarry Smith } 8773ad4599aSBarry Smith ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr); 8783ad4599aSBarry Smith if (dmhasrestrict && !mglevels[i+1]->restrct){ 8793ad4599aSBarry Smith ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr); 8803ad4599aSBarry Smith ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr); 8813ad4599aSBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 8823ad4599aSBarry Smith } 883eab52d2bSLawrence Mitchell ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr); 884eab52d2bSLawrence Mitchell if (dmhasinject && !mglevels[i+1]->inject){ 885eab52d2bSLawrence Mitchell ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr); 886eab52d2bSLawrence Mitchell ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr); 887eab52d2bSLawrence Mitchell ierr = MatDestroy(&p);CHKERRQ(ierr); 888eab52d2bSLawrence Mitchell } 88924c3aa18SBarry Smith } 8902d2b81a6SBarry Smith 891ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 8922d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 893ce4cda84SJed Brown } 894cab2e9ccSBarry Smith 895ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 8962134b1e4SBarry Smith /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 897cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 898cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 899c27ee7a3SPatrick Farrell if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) { 900c27ee7a3SPatrick Farrell ierr = KSPSetDM(mglevels[n-1]->smoothu,pc->dm);CHKERRQ(ierr); 901c27ee7a3SPatrick Farrell ierr = KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);CHKERRQ(ierr); 902c27ee7a3SPatrick Farrell } 903*41b6fd38SMatthew G. Knepley if (mglevels[n-1]->cr) { 904*41b6fd38SMatthew G. Knepley ierr = KSPSetDM(mglevels[n-1]->cr,pc->dm);CHKERRQ(ierr); 905*41b6fd38SMatthew G. Knepley ierr = KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE);CHKERRQ(ierr); 906*41b6fd38SMatthew G. Knepley } 907218a07d4SBarry Smith } 908218a07d4SBarry Smith 9092134b1e4SBarry Smith if (mg->galerkin < PC_MG_GALERKIN_NONE) { 9102134b1e4SBarry Smith Mat A,B; 9112134b1e4SBarry Smith PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 9122134b1e4SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 9132134b1e4SBarry Smith 9142134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 9152134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE; 9162134b1e4SBarry Smith if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 917f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 918b9367914SBarry 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"); 919b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 920b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 921b9367914SBarry Smith } 922b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 923b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 924b9367914SBarry Smith } 9252134b1e4SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 9262134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr); 9272134b1e4SBarry Smith } 9282134b1e4SBarry Smith if (doA) { 9292df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr); 9302134b1e4SBarry Smith } 9312134b1e4SBarry Smith if (doB) { 9322df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr); 933b9367914SBarry Smith } 9342134b1e4SBarry Smith /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 9352134b1e4SBarry Smith if (!doA && dAeqdB) { 9362134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);} 9372134b1e4SBarry Smith A = B; 9382134b1e4SBarry Smith } else if (!doA && reuse == MAT_INITIAL_MATRIX) { 9392134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr); 9402134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 941b9367914SBarry Smith } 9422134b1e4SBarry Smith if (!doB && dAeqdB) { 9432134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);} 9442134b1e4SBarry Smith B = A; 9452134b1e4SBarry Smith } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 94623ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 9472134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 9487d537d92SJed Brown } 9492134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 9502134b1e4SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr); 9512134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr); 9522134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr); 9532134b1e4SBarry Smith } 9542134b1e4SBarry Smith dA = A; 955f3fbd535SBarry Smith dB = B; 956f3fbd535SBarry Smith } 957f3fbd535SBarry Smith } 958f3b08a26SMatthew G. Knepley 959f3b08a26SMatthew G. Knepley 960f3b08a26SMatthew G. Knepley /* Adapt interpolation matrices */ 961f3b08a26SMatthew G. Knepley if (mg->adaptInterpolation) { 962f3b08a26SMatthew G. Knepley mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */ 963f3b08a26SMatthew G. Knepley for (i = 0; i < n; ++i) { 964f3b08a26SMatthew G. Knepley ierr = PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);CHKERRQ(ierr); 965f3b08a26SMatthew 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);} 966f3b08a26SMatthew G. Knepley } 967f3b08a26SMatthew G. Knepley for (i = n-2; i > -1; --i) { 968f3b08a26SMatthew G. Knepley ierr = PCMGRecomputeLevelOperators_Internal(pc, i);CHKERRQ(ierr); 969f3b08a26SMatthew G. Knepley } 970f3b08a26SMatthew G. Knepley } 971f3b08a26SMatthew G. Knepley 9722134b1e4SBarry Smith if (needRestricts && pc->dm) { 973caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 974ccceb50cSJed Brown DM dmfine,dmcoarse; 975caa4e7f2SJed Brown Mat Restrict,Inject; 976caa4e7f2SJed Brown Vec rscale; 977ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 978ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 979caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 980caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 981eab52d2bSLawrence Mitchell ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr); 982caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 983caa4e7f2SJed Brown } 984caa4e7f2SJed Brown } 985f3fbd535SBarry Smith 986f3fbd535SBarry Smith if (!pc->setupcalled) { 987f3fbd535SBarry Smith for (i=0; i<n; i++) { 988f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 989f3fbd535SBarry Smith } 990f3fbd535SBarry Smith for (i=1; i<n; i++) { 991f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 992f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 993f3fbd535SBarry Smith } 994*41b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 995*41b6fd38SMatthew G. Knepley ierr = KSPSetFromOptions(mglevels[i]->cr);CHKERRQ(ierr); 996*41b6fd38SMatthew G. Knepley } 997f3fbd535SBarry Smith } 9983ad4599aSBarry Smith /* insure that if either interpolation or restriction is set the other other one is set */ 999f3fbd535SBarry Smith for (i=1; i<n; i++) { 10003ad4599aSBarry Smith ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr); 10013ad4599aSBarry Smith ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr); 1002f3fbd535SBarry Smith } 1003f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 1004f3fbd535SBarry Smith if (!mglevels[i]->b) { 1005f3fbd535SBarry Smith Vec *vec; 10062a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 1007f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 10086bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 1009f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 1010f3fbd535SBarry Smith } 1011f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 1012f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 1013f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 10146bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 1015f3fbd535SBarry Smith } 1016f3fbd535SBarry Smith if (!mglevels[i]->x) { 1017f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 1018f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 10196bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 1020f3fbd535SBarry Smith } 1021*41b6fd38SMatthew G. Knepley if (doCR) { 1022*41b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crx);CHKERRQ(ierr); 1023*41b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crb);CHKERRQ(ierr); 1024*41b6fd38SMatthew G. Knepley ierr = VecZeroEntries(mglevels[i]->crb);CHKERRQ(ierr); 1025*41b6fd38SMatthew G. Knepley } 1026f3fbd535SBarry Smith } 1027f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 1028f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 1029f3fbd535SBarry Smith Vec *vec; 10302a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 1031f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 10326bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 1033f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 1034f3fbd535SBarry Smith } 1035*41b6fd38SMatthew G. Knepley if (doCR) { 1036*41b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx);CHKERRQ(ierr); 1037*41b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb);CHKERRQ(ierr); 1038*41b6fd38SMatthew G. Knepley ierr = VecZeroEntries(mglevels[n-1]->crb);CHKERRQ(ierr); 1039*41b6fd38SMatthew G. Knepley } 1040f3fbd535SBarry Smith } 1041f3fbd535SBarry Smith 104298e04984SBarry Smith if (pc->dm) { 104398e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 104498e04984SBarry Smith for (i=0; i<n-1; i++) { 104598e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 104698e04984SBarry Smith } 104798e04984SBarry Smith } 1048f3fbd535SBarry Smith 1049f3fbd535SBarry Smith for (i=1; i<n; i++) { 10502a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 1051f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 1052f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 1053f3fbd535SBarry Smith } 1054*41b6fd38SMatthew G. Knepley if (mglevels[i]->cr) {ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr);} 105563e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1056f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 1057c0decd05SBarry Smith if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1058899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1059899639b0SHong Zhang } 106063e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1061d42688cbSBarry Smith if (!mglevels[i]->residual) { 1062d42688cbSBarry Smith Mat mat; 106313842ffbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 106454b2cd4bSJed Brown ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 1065d42688cbSBarry Smith } 1066f3fbd535SBarry Smith } 1067f3fbd535SBarry Smith for (i=1; i<n; i++) { 1068f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 1069f3fbd535SBarry Smith Mat downmat,downpmat; 1070f3fbd535SBarry Smith 1071f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 10720298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 1073f3fbd535SBarry Smith if (!opsset) { 107423ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 107523ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 1076f3fbd535SBarry Smith } 1077f3fbd535SBarry Smith 1078f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 107963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1080f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 1081c0decd05SBarry Smith if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) { 1082899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1083899639b0SHong Zhang } 108463e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1085f3fbd535SBarry Smith } 1086*41b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 1087*41b6fd38SMatthew G. Knepley Mat downmat,downpmat; 1088*41b6fd38SMatthew G. Knepley 1089*41b6fd38SMatthew G. Knepley /* check if operators have been set for up, if not use down operators to set them */ 1090*41b6fd38SMatthew G. Knepley ierr = KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL);CHKERRQ(ierr); 1091*41b6fd38SMatthew G. Knepley if (!opsset) { 1092*41b6fd38SMatthew G. Knepley ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 1093*41b6fd38SMatthew G. Knepley ierr = KSPSetOperators(mglevels[i]->cr,downmat,downpmat);CHKERRQ(ierr); 1094*41b6fd38SMatthew G. Knepley } 1095*41b6fd38SMatthew G. Knepley 1096*41b6fd38SMatthew G. Knepley ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr); 1097*41b6fd38SMatthew G. Knepley if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1098*41b6fd38SMatthew G. Knepley ierr = KSPSetUp(mglevels[i]->cr);CHKERRQ(ierr); 1099*41b6fd38SMatthew G. Knepley if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) { 1100*41b6fd38SMatthew G. Knepley pc->failedreason = PC_SUBPC_ERROR; 1101*41b6fd38SMatthew G. Knepley } 1102*41b6fd38SMatthew G. Knepley if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1103*41b6fd38SMatthew G. Knepley } 1104f3fbd535SBarry Smith } 1105f3fbd535SBarry Smith 110663e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1107f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 1108c0decd05SBarry Smith if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1109899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1110899639b0SHong Zhang } 111163e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1112f3fbd535SBarry Smith 1113f3fbd535SBarry Smith /* 1114f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 1115e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 1116f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 1117f3fbd535SBarry Smith 1118f3fbd535SBarry Smith Only support one or the other at the same time. 1119f3fbd535SBarry Smith */ 1120f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 1121c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 1122ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 1123f3fbd535SBarry Smith dump = PETSC_FALSE; 1124f3fbd535SBarry Smith #endif 1125c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 1126ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 1127f3fbd535SBarry Smith 1128f3fbd535SBarry Smith if (viewer) { 1129f3fbd535SBarry Smith for (i=1; i<n; i++) { 1130f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 1131f3fbd535SBarry Smith } 1132f3fbd535SBarry Smith for (i=0; i<n; i++) { 1133f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 1134f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 1135f3fbd535SBarry Smith } 1136f3fbd535SBarry Smith } 1137f3fbd535SBarry Smith PetscFunctionReturn(0); 1138f3fbd535SBarry Smith } 1139f3fbd535SBarry Smith 1140f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 1141f3fbd535SBarry Smith 114297d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 114397d33e41SMatthew G. Knepley { 114497d33e41SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 114597d33e41SMatthew G. Knepley 114697d33e41SMatthew G. Knepley PetscFunctionBegin; 114797d33e41SMatthew G. Knepley *levels = mg->nlevels; 114897d33e41SMatthew G. Knepley PetscFunctionReturn(0); 114997d33e41SMatthew G. Knepley } 115097d33e41SMatthew G. Knepley 11514b9ad928SBarry Smith /*@ 115297177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 11534b9ad928SBarry Smith 11544b9ad928SBarry Smith Not Collective 11554b9ad928SBarry Smith 11564b9ad928SBarry Smith Input Parameter: 11574b9ad928SBarry Smith . pc - the preconditioner context 11584b9ad928SBarry Smith 11594b9ad928SBarry Smith Output parameter: 11604b9ad928SBarry Smith . levels - the number of levels 11614b9ad928SBarry Smith 11624b9ad928SBarry Smith Level: advanced 11634b9ad928SBarry Smith 116497177400SBarry Smith .seealso: PCMGSetLevels() 11654b9ad928SBarry Smith @*/ 11667087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 11674b9ad928SBarry Smith { 1168e952c529SMatthew G. Knepley PetscErrorCode ierr; 11694b9ad928SBarry Smith 11704b9ad928SBarry Smith PetscFunctionBegin; 11710700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11724482741eSBarry Smith PetscValidIntPointer(levels,2); 117397d33e41SMatthew G. Knepley *levels = 0; 117437b5128cSMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr); 11754b9ad928SBarry Smith PetscFunctionReturn(0); 11764b9ad928SBarry Smith } 11774b9ad928SBarry Smith 11784b9ad928SBarry Smith /*@ 117997177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 11804b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 11814b9ad928SBarry Smith 1182ad4df100SBarry Smith Logically Collective on PC 11834b9ad928SBarry Smith 11844b9ad928SBarry Smith Input Parameters: 11854b9ad928SBarry Smith + pc - the preconditioner context 11869dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 11879dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 11884b9ad928SBarry Smith 11894b9ad928SBarry Smith Options Database Key: 11904b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 11914b9ad928SBarry Smith additive, full, kaskade 11924b9ad928SBarry Smith 11934b9ad928SBarry Smith Level: advanced 11944b9ad928SBarry Smith 119597177400SBarry Smith .seealso: PCMGSetLevels() 11964b9ad928SBarry Smith @*/ 11977087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 11984b9ad928SBarry Smith { 1199f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 12004b9ad928SBarry Smith 12014b9ad928SBarry Smith PetscFunctionBegin; 12020700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1203c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 120431567311SBarry Smith mg->am = form; 12059dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 1206c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 1207c60c7ad4SBarry Smith PetscFunctionReturn(0); 1208c60c7ad4SBarry Smith } 1209c60c7ad4SBarry Smith 1210c60c7ad4SBarry Smith /*@ 1211c60c7ad4SBarry Smith PCMGGetType - Determines the form of multigrid to use: 1212c60c7ad4SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 1213c60c7ad4SBarry Smith 1214c60c7ad4SBarry Smith Logically Collective on PC 1215c60c7ad4SBarry Smith 1216c60c7ad4SBarry Smith Input Parameter: 1217c60c7ad4SBarry Smith . pc - the preconditioner context 1218c60c7ad4SBarry Smith 1219c60c7ad4SBarry Smith Output Parameter: 1220c60c7ad4SBarry Smith . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 1221c60c7ad4SBarry Smith 1222c60c7ad4SBarry Smith 1223c60c7ad4SBarry Smith Level: advanced 1224c60c7ad4SBarry Smith 1225c60c7ad4SBarry Smith .seealso: PCMGSetLevels() 1226c60c7ad4SBarry Smith @*/ 1227c60c7ad4SBarry Smith PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 1228c60c7ad4SBarry Smith { 1229c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1230c60c7ad4SBarry Smith 1231c60c7ad4SBarry Smith PetscFunctionBegin; 1232c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1233c60c7ad4SBarry Smith *type = mg->am; 12344b9ad928SBarry Smith PetscFunctionReturn(0); 12354b9ad928SBarry Smith } 12364b9ad928SBarry Smith 12374b9ad928SBarry Smith /*@ 12380d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 12394b9ad928SBarry Smith complicated cycling. 12404b9ad928SBarry Smith 1241ad4df100SBarry Smith Logically Collective on PC 12424b9ad928SBarry Smith 12434b9ad928SBarry Smith Input Parameters: 1244c2be2410SBarry Smith + pc - the multigrid context 1245c1cbb1deSBarry Smith - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 12464b9ad928SBarry Smith 12474b9ad928SBarry Smith Options Database Key: 1248c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired 12494b9ad928SBarry Smith 12504b9ad928SBarry Smith Level: advanced 12514b9ad928SBarry Smith 12520d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 12534b9ad928SBarry Smith @*/ 12547087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 12554b9ad928SBarry Smith { 1256f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1257f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 125879416396SBarry Smith PetscInt i,levels; 12594b9ad928SBarry Smith 12604b9ad928SBarry Smith PetscFunctionBegin; 12610700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 126269fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc,n,2); 12631a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1264f3fbd535SBarry Smith levels = mglevels[0]->levels; 12652fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 12664b9ad928SBarry Smith PetscFunctionReturn(0); 12674b9ad928SBarry Smith } 12684b9ad928SBarry Smith 12698cc2d5dfSBarry Smith /*@ 12708cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 12718cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 12728cc2d5dfSBarry Smith 1273ad4df100SBarry Smith Logically Collective on PC 12748cc2d5dfSBarry Smith 12758cc2d5dfSBarry Smith Input Parameters: 12768cc2d5dfSBarry Smith + pc - the multigrid context 12778cc2d5dfSBarry Smith - n - number of cycles (default is 1) 12788cc2d5dfSBarry Smith 12798cc2d5dfSBarry Smith Options Database Key: 1280e1bc860dSBarry Smith . -pc_mg_multiplicative_cycles n 12818cc2d5dfSBarry Smith 12828cc2d5dfSBarry Smith Level: advanced 12838cc2d5dfSBarry Smith 128495452b02SPatrick Sanan Notes: 128595452b02SPatrick Sanan This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 12868cc2d5dfSBarry Smith 12878cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 12888cc2d5dfSBarry Smith @*/ 12897087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 12908cc2d5dfSBarry Smith { 1291f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 12928cc2d5dfSBarry Smith 12938cc2d5dfSBarry Smith PetscFunctionBegin; 12940700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1295c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 12962a21e185SBarry Smith mg->cyclesperpcapply = n; 12978cc2d5dfSBarry Smith PetscFunctionReturn(0); 12988cc2d5dfSBarry Smith } 12998cc2d5dfSBarry Smith 13002134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1301fb15c04eSBarry Smith { 1302fb15c04eSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1303fb15c04eSBarry Smith 1304fb15c04eSBarry Smith PetscFunctionBegin; 13052134b1e4SBarry Smith mg->galerkin = use; 1306fb15c04eSBarry Smith PetscFunctionReturn(0); 1307fb15c04eSBarry Smith } 1308fb15c04eSBarry Smith 1309c2be2410SBarry Smith /*@ 131097177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 131182b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1312c2be2410SBarry Smith 1313ad4df100SBarry Smith Logically Collective on PC 1314c2be2410SBarry Smith 1315c2be2410SBarry Smith Input Parameters: 1316c91913d3SJed Brown + pc - the multigrid context 13172134b1e4SBarry Smith - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1318c2be2410SBarry Smith 1319c2be2410SBarry Smith Options Database Key: 13202134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> 1321c2be2410SBarry Smith 1322c2be2410SBarry Smith Level: intermediate 1323c2be2410SBarry Smith 132495452b02SPatrick Sanan Notes: 132595452b02SPatrick Sanan Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 13261c1aac46SBarry Smith use the PCMG construction of the coarser grids. 13271c1aac46SBarry Smith 13282134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType 13293fc8bf9cSBarry Smith 1330c2be2410SBarry Smith @*/ 13312134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1332c2be2410SBarry Smith { 1333fb15c04eSBarry Smith PetscErrorCode ierr; 1334c2be2410SBarry Smith 1335c2be2410SBarry Smith PetscFunctionBegin; 13360700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13372134b1e4SBarry Smith ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr); 1338c2be2410SBarry Smith PetscFunctionReturn(0); 1339c2be2410SBarry Smith } 1340c2be2410SBarry Smith 13413fc8bf9cSBarry Smith /*@ 13423fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 134382b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 13443fc8bf9cSBarry Smith 13453fc8bf9cSBarry Smith Not Collective 13463fc8bf9cSBarry Smith 13473fc8bf9cSBarry Smith Input Parameter: 13483fc8bf9cSBarry Smith . pc - the multigrid context 13493fc8bf9cSBarry Smith 13503fc8bf9cSBarry Smith Output Parameter: 13512134b1e4SBarry 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 13523fc8bf9cSBarry Smith 13533fc8bf9cSBarry Smith Level: intermediate 13543fc8bf9cSBarry Smith 13552134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType 13563fc8bf9cSBarry Smith 13573fc8bf9cSBarry Smith @*/ 13582134b1e4SBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 13593fc8bf9cSBarry Smith { 1360f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 13613fc8bf9cSBarry Smith 13623fc8bf9cSBarry Smith PetscFunctionBegin; 13630700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13642134b1e4SBarry Smith *galerkin = mg->galerkin; 13653fc8bf9cSBarry Smith PetscFunctionReturn(0); 13663fc8bf9cSBarry Smith } 13673fc8bf9cSBarry Smith 1368f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt) 1369f3b08a26SMatthew G. Knepley { 1370f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1371f3b08a26SMatthew G. Knepley 1372f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1373f3b08a26SMatthew G. Knepley mg->adaptInterpolation = adapt; 1374f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1375f3b08a26SMatthew G. Knepley } 1376f3b08a26SMatthew G. Knepley 1377f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt) 1378f3b08a26SMatthew G. Knepley { 1379f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1380f3b08a26SMatthew G. Knepley 1381f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1382f3b08a26SMatthew G. Knepley *adapt = mg->adaptInterpolation; 1383f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1384f3b08a26SMatthew G. Knepley } 1385f3b08a26SMatthew G. Knepley 1386*41b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr) 1387*41b6fd38SMatthew G. Knepley { 1388*41b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1389*41b6fd38SMatthew G. Knepley 1390*41b6fd38SMatthew G. Knepley PetscFunctionBegin; 1391*41b6fd38SMatthew G. Knepley mg->compatibleRelaxation = cr; 1392*41b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 1393*41b6fd38SMatthew G. Knepley } 1394*41b6fd38SMatthew G. Knepley 1395*41b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr) 1396*41b6fd38SMatthew G. Knepley { 1397*41b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1398*41b6fd38SMatthew G. Knepley 1399*41b6fd38SMatthew G. Knepley PetscFunctionBegin; 1400*41b6fd38SMatthew G. Knepley *cr = mg->compatibleRelaxation; 1401*41b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 1402*41b6fd38SMatthew G. Knepley } 1403*41b6fd38SMatthew G. Knepley 1404f3b08a26SMatthew G. Knepley /*@ 1405f3b08a26SMatthew 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. 1406f3b08a26SMatthew G. Knepley 1407f3b08a26SMatthew G. Knepley Logically Collective on PC 1408f3b08a26SMatthew G. Knepley 1409f3b08a26SMatthew G. Knepley Input Parameters: 1410f3b08a26SMatthew G. Knepley + pc - the multigrid context 1411f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator 1412f3b08a26SMatthew G. Knepley 1413f3b08a26SMatthew G. Knepley Options Database Keys: 1414f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp - Turn on adaptation 1415f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension 1416f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector 1417f3b08a26SMatthew G. Knepley 1418f3b08a26SMatthew G. Knepley Level: intermediate 1419f3b08a26SMatthew G. Knepley 1420f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1421f3b08a26SMatthew G. Knepley .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin() 1422f3b08a26SMatthew G. Knepley @*/ 1423f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt) 1424f3b08a26SMatthew G. Knepley { 1425f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1426f3b08a26SMatthew G. Knepley 1427f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1428f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1429f3b08a26SMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));CHKERRQ(ierr); 1430f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1431f3b08a26SMatthew G. Knepley } 1432f3b08a26SMatthew G. Knepley 1433f3b08a26SMatthew G. Knepley /*@ 1434f3b08a26SMatthew 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. 1435f3b08a26SMatthew G. Knepley 1436f3b08a26SMatthew G. Knepley Not collective 1437f3b08a26SMatthew G. Knepley 1438f3b08a26SMatthew G. Knepley Input Parameter: 1439f3b08a26SMatthew G. Knepley . pc - the multigrid context 1440f3b08a26SMatthew G. Knepley 1441f3b08a26SMatthew G. Knepley Output Parameter: 1442f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator 1443f3b08a26SMatthew G. Knepley 1444f3b08a26SMatthew G. Knepley Level: intermediate 1445f3b08a26SMatthew G. Knepley 1446f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1447f3b08a26SMatthew G. Knepley .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin() 1448f3b08a26SMatthew G. Knepley @*/ 1449f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt) 1450f3b08a26SMatthew G. Knepley { 1451f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1452f3b08a26SMatthew G. Knepley 1453f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1454f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1455f3b08a26SMatthew G. Knepley PetscValidPointer(adapt, 2); 1456f3b08a26SMatthew G. Knepley ierr = PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));CHKERRQ(ierr); 1457f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1458f3b08a26SMatthew G. Knepley } 1459f3b08a26SMatthew G. Knepley 14604b9ad928SBarry Smith /*@ 1461*41b6fd38SMatthew G. Knepley PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation. 1462*41b6fd38SMatthew G. Knepley 1463*41b6fd38SMatthew G. Knepley Logically Collective on PC 1464*41b6fd38SMatthew G. Knepley 1465*41b6fd38SMatthew G. Knepley Input Parameters: 1466*41b6fd38SMatthew G. Knepley + pc - the multigrid context 1467*41b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation 1468*41b6fd38SMatthew G. Knepley 1469*41b6fd38SMatthew G. Knepley Options Database Keys: 1470*41b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation 1471*41b6fd38SMatthew G. Knepley 1472*41b6fd38SMatthew G. Knepley Level: intermediate 1473*41b6fd38SMatthew G. Knepley 1474*41b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin 1475*41b6fd38SMatthew G. Knepley .seealso: PCMGGetAdaptCR(), PCMGSetAdaptInterpolation(), PCMGSetGalerkin() 1476*41b6fd38SMatthew G. Knepley @*/ 1477*41b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr) 1478*41b6fd38SMatthew G. Knepley { 1479*41b6fd38SMatthew G. Knepley PetscErrorCode ierr; 1480*41b6fd38SMatthew G. Knepley 1481*41b6fd38SMatthew G. Knepley PetscFunctionBegin; 1482*41b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1483*41b6fd38SMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr));CHKERRQ(ierr); 1484*41b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 1485*41b6fd38SMatthew G. Knepley } 1486*41b6fd38SMatthew G. Knepley 1487*41b6fd38SMatthew G. Knepley /*@ 1488*41b6fd38SMatthew G. Knepley PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation. 1489*41b6fd38SMatthew G. Knepley 1490*41b6fd38SMatthew G. Knepley Not collective 1491*41b6fd38SMatthew G. Knepley 1492*41b6fd38SMatthew G. Knepley Input Parameter: 1493*41b6fd38SMatthew G. Knepley . pc - the multigrid context 1494*41b6fd38SMatthew G. Knepley 1495*41b6fd38SMatthew G. Knepley Output Parameter: 1496*41b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion 1497*41b6fd38SMatthew G. Knepley 1498*41b6fd38SMatthew G. Knepley Level: intermediate 1499*41b6fd38SMatthew G. Knepley 1500*41b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin 1501*41b6fd38SMatthew G. Knepley .seealso: PCMGSetAdaptCR(), PCMGGetAdaptInterpolation(), PCMGSetGalerkin() 1502*41b6fd38SMatthew G. Knepley @*/ 1503*41b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr) 1504*41b6fd38SMatthew G. Knepley { 1505*41b6fd38SMatthew G. Knepley PetscErrorCode ierr; 1506*41b6fd38SMatthew G. Knepley 1507*41b6fd38SMatthew G. Knepley PetscFunctionBegin; 1508*41b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1509*41b6fd38SMatthew G. Knepley PetscValidPointer(cr, 2); 1510*41b6fd38SMatthew G. Knepley ierr = PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr));CHKERRQ(ierr); 1511*41b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 1512*41b6fd38SMatthew G. Knepley } 1513*41b6fd38SMatthew G. Knepley 1514*41b6fd38SMatthew G. Knepley /*@ 151506792cafSBarry Smith PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1516710315b6SLawrence Mitchell on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1517710315b6SLawrence Mitchell pre- and post-smoothing steps. 151806792cafSBarry Smith 151906792cafSBarry Smith Logically Collective on PC 152006792cafSBarry Smith 152106792cafSBarry Smith Input Parameters: 152206792cafSBarry Smith + mg - the multigrid context 152306792cafSBarry Smith - n - the number of smoothing steps 152406792cafSBarry Smith 152506792cafSBarry Smith Options Database Key: 1526a2b725a8SWilliam Gropp . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 152706792cafSBarry Smith 152806792cafSBarry Smith Level: advanced 152906792cafSBarry Smith 153095452b02SPatrick Sanan Notes: 153195452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 153206792cafSBarry Smith there is no separate smooth up on the coarsest grid. 153306792cafSBarry Smith 1534710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp() 153506792cafSBarry Smith @*/ 153606792cafSBarry Smith PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 153706792cafSBarry Smith { 153806792cafSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 153906792cafSBarry Smith PC_MG_Levels **mglevels = mg->levels; 154006792cafSBarry Smith PetscErrorCode ierr; 154106792cafSBarry Smith PetscInt i,levels; 154206792cafSBarry Smith 154306792cafSBarry Smith PetscFunctionBegin; 154406792cafSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 154506792cafSBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 15461a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 154706792cafSBarry Smith levels = mglevels[0]->levels; 154806792cafSBarry Smith 154906792cafSBarry Smith for (i=1; i<levels; i++) { 155006792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 155106792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 155206792cafSBarry Smith mg->default_smoothu = n; 155306792cafSBarry Smith mg->default_smoothd = n; 155406792cafSBarry Smith } 155506792cafSBarry Smith PetscFunctionReturn(0); 155606792cafSBarry Smith } 155706792cafSBarry Smith 1558f442ab6aSBarry Smith /*@ 15598e5aa403SBarry Smith PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels 1560710315b6SLawrence Mitchell and adds the suffix _up to the options name 1561f442ab6aSBarry Smith 1562f442ab6aSBarry Smith Logically Collective on PC 1563f442ab6aSBarry Smith 1564f442ab6aSBarry Smith Input Parameters: 1565f442ab6aSBarry Smith . pc - the preconditioner context 1566f442ab6aSBarry Smith 1567f442ab6aSBarry Smith Options Database Key: 1568f442ab6aSBarry Smith . -pc_mg_distinct_smoothup 1569f442ab6aSBarry Smith 1570f442ab6aSBarry Smith Level: advanced 1571f442ab6aSBarry Smith 157295452b02SPatrick Sanan Notes: 157395452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 1574f442ab6aSBarry Smith there is no separate smooth up on the coarsest grid. 1575f442ab6aSBarry Smith 1576710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth() 1577f442ab6aSBarry Smith @*/ 1578f442ab6aSBarry Smith PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1579f442ab6aSBarry Smith { 1580f442ab6aSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1581f442ab6aSBarry Smith PC_MG_Levels **mglevels = mg->levels; 1582f442ab6aSBarry Smith PetscErrorCode ierr; 1583f442ab6aSBarry Smith PetscInt i,levels; 1584f442ab6aSBarry Smith KSP subksp; 1585f442ab6aSBarry Smith 1586f442ab6aSBarry Smith PetscFunctionBegin; 1587f442ab6aSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1588f442ab6aSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1589f442ab6aSBarry Smith levels = mglevels[0]->levels; 1590f442ab6aSBarry Smith 1591f442ab6aSBarry Smith for (i=1; i<levels; i++) { 1592710315b6SLawrence Mitchell const char *prefix = NULL; 1593f442ab6aSBarry Smith /* make sure smoother up and down are different */ 1594f442ab6aSBarry Smith ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr); 1595710315b6SLawrence Mitchell ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr); 1596710315b6SLawrence Mitchell ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr); 1597f442ab6aSBarry Smith ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr); 1598f442ab6aSBarry Smith } 1599f442ab6aSBarry Smith PetscFunctionReturn(0); 1600f442ab6aSBarry Smith } 1601f442ab6aSBarry Smith 160207a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 160307a4832bSFande Kong PetscErrorCode PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[]) 160407a4832bSFande Kong { 160507a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 160607a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 160707a4832bSFande Kong Mat *mat; 160807a4832bSFande Kong PetscInt l; 160907a4832bSFande Kong PetscErrorCode ierr; 161007a4832bSFande Kong 161107a4832bSFande Kong PetscFunctionBegin; 161207a4832bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 161307a4832bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 161407a4832bSFande Kong for (l=1; l< mg->nlevels; l++) { 161507a4832bSFande Kong mat[l-1] = mglevels[l]->interpolate; 161607a4832bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr); 161707a4832bSFande Kong } 161807a4832bSFande Kong *num_levels = mg->nlevels; 161907a4832bSFande Kong *interpolations = mat; 162007a4832bSFande Kong PetscFunctionReturn(0); 162107a4832bSFande Kong } 162207a4832bSFande Kong 162307a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 162407a4832bSFande Kong PetscErrorCode PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[]) 162507a4832bSFande Kong { 162607a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 162707a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 162807a4832bSFande Kong PetscInt l; 162907a4832bSFande Kong Mat *mat; 163007a4832bSFande Kong PetscErrorCode ierr; 163107a4832bSFande Kong 163207a4832bSFande Kong PetscFunctionBegin; 163307a4832bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 163407a4832bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 163507a4832bSFande Kong for (l=0; l<mg->nlevels-1; l++) { 163607a4832bSFande Kong ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr); 163707a4832bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr); 163807a4832bSFande Kong } 163907a4832bSFande Kong *num_levels = mg->nlevels; 164007a4832bSFande Kong *coarseOperators = mat; 164107a4832bSFande Kong PetscFunctionReturn(0); 164207a4832bSFande Kong } 164307a4832bSFande Kong 1644f3b08a26SMatthew G. Knepley /*@C 1645f3b08a26SMatthew G. Knepley PCMGRegisterCoarseSpaceConstructor - Adds a method to the PCMG package for coarse space construction. 1646f3b08a26SMatthew G. Knepley 1647f3b08a26SMatthew G. Knepley Not collective 1648f3b08a26SMatthew G. Knepley 1649f3b08a26SMatthew G. Knepley Input Parameters: 1650f3b08a26SMatthew G. Knepley + name - name of the constructor 1651f3b08a26SMatthew G. Knepley - function - constructor routine 1652f3b08a26SMatthew G. Knepley 1653f3b08a26SMatthew G. Knepley Notes: 1654f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1655f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp) 1656f3b08a26SMatthew G. Knepley $ pc - The PC object 1657f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1658f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1659f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1660f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1661f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1662f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1663f3b08a26SMatthew G. Knepley 1664f3b08a26SMatthew G. Knepley Level: advanced 1665f3b08a26SMatthew G. Knepley 1666f3b08a26SMatthew G. Knepley .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister() 1667f3b08a26SMatthew G. Knepley @*/ 1668f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)) 1669f3b08a26SMatthew G. Knepley { 1670f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1671f3b08a26SMatthew G. Knepley 1672f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1673f3b08a26SMatthew G. Knepley ierr = PCInitializePackage();CHKERRQ(ierr); 1674f3b08a26SMatthew G. Knepley ierr = PetscFunctionListAdd(&PCMGCoarseList,name,function);CHKERRQ(ierr); 1675f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1676f3b08a26SMatthew G. Knepley } 1677f3b08a26SMatthew G. Knepley 1678f3b08a26SMatthew G. Knepley /*@C 1679f3b08a26SMatthew G. Knepley PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method. 1680f3b08a26SMatthew G. Knepley 1681f3b08a26SMatthew G. Knepley Not collective 1682f3b08a26SMatthew G. Knepley 1683f3b08a26SMatthew G. Knepley Input Parameter: 1684f3b08a26SMatthew G. Knepley . name - name of the constructor 1685f3b08a26SMatthew G. Knepley 1686f3b08a26SMatthew G. Knepley Output Parameter: 1687f3b08a26SMatthew G. Knepley . function - constructor routine 1688f3b08a26SMatthew G. Knepley 1689f3b08a26SMatthew G. Knepley Notes: 1690f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1691f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp) 1692f3b08a26SMatthew G. Knepley $ pc - The PC object 1693f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1694f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1695f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1696f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1697f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1698f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1699f3b08a26SMatthew G. Knepley 1700f3b08a26SMatthew G. Knepley Level: advanced 1701f3b08a26SMatthew G. Knepley 1702f3b08a26SMatthew G. Knepley .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister() 1703f3b08a26SMatthew G. Knepley @*/ 1704f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)) 1705f3b08a26SMatthew G. Knepley { 1706f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1707f3b08a26SMatthew G. Knepley 1708f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1709f3b08a26SMatthew G. Knepley ierr = PetscFunctionListFind(PCMGCoarseList,name,function);CHKERRQ(ierr); 1710f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1711f3b08a26SMatthew G. Knepley } 1712f3b08a26SMatthew G. Knepley 17134b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 17144b9ad928SBarry Smith 17153b09bd56SBarry Smith /*MC 1716ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 17173b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 17183b09bd56SBarry Smith 17193b09bd56SBarry Smith Options Database Keys: 17203b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1721391689abSStefano Zampini . -pc_mg_cycle_type <v,w> - provide the cycle desired 17228c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 17233b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 1724710315b6SLawrence Mitchell . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 17252134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 17268cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 17278cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1728e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 17298cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 17308cf6037eSBarry Smith to the binary output file called binaryoutput 17313b09bd56SBarry Smith 173295452b02SPatrick Sanan Notes: 17330b3c753eSRichard 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 17343b09bd56SBarry Smith 17358cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 17368cf6037eSBarry Smith 173723067569SBarry Smith When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 173823067569SBarry Smith is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 173923067569SBarry Smith (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 174023067569SBarry Smith residual is computed at the end of each cycle. 174123067569SBarry Smith 17423b09bd56SBarry Smith Level: intermediate 17433b09bd56SBarry Smith 17448cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1745710315b6SLawrence Mitchell PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), 1746710315b6SLawrence Mitchell PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 174797177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 17480d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 17493b09bd56SBarry Smith M*/ 17503b09bd56SBarry Smith 17518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 17524b9ad928SBarry Smith { 1753f3fbd535SBarry Smith PC_MG *mg; 1754f3fbd535SBarry Smith PetscErrorCode ierr; 1755f3fbd535SBarry Smith 17564b9ad928SBarry Smith PetscFunctionBegin; 1757b00a9115SJed Brown ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1758f3fbd535SBarry Smith pc->data = (void*)mg; 1759f3fbd535SBarry Smith mg->nlevels = -1; 176010eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 17612134b1e4SBarry Smith mg->galerkin = PC_MG_GALERKIN_NONE; 1762f3b08a26SMatthew G. Knepley mg->adaptInterpolation = PETSC_FALSE; 1763f3b08a26SMatthew G. Knepley mg->Nc = -1; 1764f3b08a26SMatthew G. Knepley mg->eigenvalue = -1; 1765f3fbd535SBarry Smith 176637a44384SMark Adams pc->useAmat = PETSC_TRUE; 176737a44384SMark Adams 17684b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 17694b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1770a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 17714b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 17724b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 17734b9ad928SBarry Smith pc->ops->view = PCView_MG; 1774fb15c04eSBarry Smith 1775f3b08a26SMatthew G. Knepley ierr = PetscObjectComposedDataRegister(&mg->eigenvalue);CHKERRQ(ierr); 1776fb15c04eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 177797d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr); 177897d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr); 1779fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);CHKERRQ(ierr); 1780fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);CHKERRQ(ierr); 1781f3b08a26SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);CHKERRQ(ierr); 1782f3b08a26SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);CHKERRQ(ierr); 1783*41b6fd38SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG);CHKERRQ(ierr); 1784*41b6fd38SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG);CHKERRQ(ierr); 17854b9ad928SBarry Smith PetscFunctionReturn(0); 17864b9ad928SBarry Smith } 1787