1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Defines the multigrid preconditioner interface. 44b9ad928SBarry Smith */ 5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> /*I "petscksp.h" I*/ 641b6fd38SMatthew 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 1530b0564aSStefano Zampini PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PetscBool transpose,PetscBool matapp,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);} 24fcb023d4SJed Brown if (!transpose) { 2530b0564aSStefano Zampini if (matapp) { 2630b0564aSStefano Zampini ierr = KSPMatSolve(mglevels->smoothd,mglevels->B,mglevels->X);CHKERRQ(ierr); /* pre-smooth */ 2730b0564aSStefano Zampini ierr = KSPCheckSolve(mglevels->smoothd,pc,NULL);CHKERRQ(ierr); 2830b0564aSStefano Zampini } else { 2931567311SBarry Smith ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ 30c0decd05SBarry Smith ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr); 3130b0564aSStefano Zampini } 32fcb023d4SJed Brown } else { 3330b0564aSStefano Zampini if (matapp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 34fcb023d4SJed Brown ierr = KSPSolveTranspose(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* transpose of post-smooth */ 35fcb023d4SJed Brown ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr); 36fcb023d4SJed Brown } 3763e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 3831567311SBarry Smith if (mglevels->level) { /* not the coarsest grid */ 3963e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 4030b0564aSStefano Zampini if (matapp && !mglevels->R) { 4130b0564aSStefano Zampini ierr = MatDuplicate(mglevels->B,MAT_DO_NOT_COPY_VALUES,&mglevels->R);CHKERRQ(ierr); 4230b0564aSStefano Zampini } 43fcb023d4SJed Brown if (!transpose) { 4430b0564aSStefano Zampini if (matapp) { ierr = (*mglevels->matresidual)(mglevels->A,mglevels->B,mglevels->X,mglevels->R);CHKERRQ(ierr); } 4530b0564aSStefano Zampini else { ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); } 46fcb023d4SJed Brown } else { 4730b0564aSStefano Zampini if (matapp) { ierr = (*mglevels->matresidualtranspose)(mglevels->A,mglevels->B,mglevels->X,mglevels->R);CHKERRQ(ierr); } 4830b0564aSStefano Zampini else { ierr = (*mglevels->residualtranspose)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); } 49fcb023d4SJed Brown } 5063e6d426SJed Brown if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} 514b9ad928SBarry Smith 524b9ad928SBarry Smith /* if on finest level and have convergence criteria set */ 5331567311SBarry Smith if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { 544b9ad928SBarry Smith PetscReal rnorm; 5531567311SBarry Smith ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); 564b9ad928SBarry Smith if (rnorm <= mg->ttol) { 5770441072SBarry Smith if (rnorm < mg->abstol) { 584d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_ATOL; 5957622a8eSBarry 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); 604b9ad928SBarry Smith } else { 614d0a8057SBarry Smith *reason = PCRICHARDSON_CONVERGED_RTOL; 6257622a8eSBarry 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); 634b9ad928SBarry Smith } 644b9ad928SBarry Smith PetscFunctionReturn(0); 654b9ad928SBarry Smith } 664b9ad928SBarry Smith } 674b9ad928SBarry Smith 6831567311SBarry Smith mgc = *(mglevelsin - 1); 6963e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 70fcb023d4SJed Brown if (!transpose) { 7130b0564aSStefano Zampini if (matapp) { ierr = MatMatRestrict(mglevels->restrct,mglevels->R,&mgc->B);CHKERRQ(ierr); } 7230b0564aSStefano Zampini else { ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); } 73fcb023d4SJed Brown } else { 7430b0564aSStefano Zampini if (matapp) { ierr = MatMatRestrict(mglevels->interpolate,mglevels->R,&mgc->B);CHKERRQ(ierr); } 7530b0564aSStefano Zampini else { ierr = MatRestrict(mglevels->interpolate,mglevels->r,mgc->b);CHKERRQ(ierr); } 76fcb023d4SJed Brown } 7763e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 7830b0564aSStefano Zampini if (matapp) { 7930b0564aSStefano Zampini if (!mgc->X) { 8030b0564aSStefano Zampini ierr = MatDuplicate(mgc->B,MAT_DO_NOT_COPY_VALUES,&mgc->X);CHKERRQ(ierr); 8130b0564aSStefano Zampini } else { 8230b0564aSStefano Zampini ierr = MatZeroEntries(mgc->X);CHKERRQ(ierr); 8330b0564aSStefano Zampini } 8430b0564aSStefano Zampini } else { 8530b0564aSStefano Zampini ierr = VecZeroEntries(mgc->x);CHKERRQ(ierr); 8630b0564aSStefano Zampini } 874b9ad928SBarry Smith while (cycles--) { 8830b0564aSStefano Zampini ierr = PCMGMCycle_Private(pc,mglevelsin-1,transpose,matapp,reason);CHKERRQ(ierr); 894b9ad928SBarry Smith } 9063e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 91fcb023d4SJed Brown if (!transpose) { 9230b0564aSStefano Zampini if (matapp) { ierr = MatMatInterpolateAdd(mglevels->interpolate,mgc->X,mglevels->X,&mglevels->X);CHKERRQ(ierr); } 9330b0564aSStefano Zampini else { ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); } 94fcb023d4SJed Brown } else { 95fcb023d4SJed Brown ierr = MatInterpolateAdd(mglevels->restrct,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); 96fcb023d4SJed Brown } 9763e6d426SJed Brown if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 9863e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 99fcb023d4SJed Brown if (!transpose) { 10030b0564aSStefano Zampini if (matapp) { 10130b0564aSStefano Zampini ierr = KSPMatSolve(mglevels->smoothu,mglevels->B,mglevels->X);CHKERRQ(ierr); /* post smooth */ 10230b0564aSStefano Zampini ierr = KSPCheckSolve(mglevels->smoothu,pc,NULL);CHKERRQ(ierr); 10330b0564aSStefano Zampini } else { 10431567311SBarry Smith ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 105c0decd05SBarry Smith ierr = KSPCheckSolve(mglevels->smoothu,pc,mglevels->x);CHKERRQ(ierr); 10630b0564aSStefano Zampini } 107fcb023d4SJed Brown } else { 10830b0564aSStefano Zampini if (matapp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 109fcb023d4SJed Brown ierr = KSPSolveTranspose(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ 110fcb023d4SJed Brown ierr = KSPCheckSolve(mglevels->smoothd,pc,mglevels->x);CHKERRQ(ierr); 111fcb023d4SJed Brown } 11241b6fd38SMatthew G. Knepley if (mglevels->cr) { 11330b0564aSStefano Zampini if (matapp) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 11441b6fd38SMatthew G. Knepley /* TODO Turn on copy and turn off noisy if we have an exact solution 11541b6fd38SMatthew G. Knepley ierr = VecCopy(mglevels->x, mglevels->crx);CHKERRQ(ierr); 11641b6fd38SMatthew G. Knepley ierr = VecCopy(mglevels->b, mglevels->crb);CHKERRQ(ierr); */ 11741b6fd38SMatthew G. Knepley ierr = KSPSetNoisy_Private(mglevels->crx);CHKERRQ(ierr); 11841b6fd38SMatthew G. Knepley ierr = KSPSolve(mglevels->cr,mglevels->crb,mglevels->crx);CHKERRQ(ierr); /* compatible relaxation */ 11941b6fd38SMatthew G. Knepley ierr = KSPCheckSolve(mglevels->cr,pc,mglevels->crx);CHKERRQ(ierr); 12041b6fd38SMatthew G. Knepley } 12163e6d426SJed Brown if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 1224b9ad928SBarry Smith } 1234b9ad928SBarry Smith PetscFunctionReturn(0); 1244b9ad928SBarry Smith } 1254b9ad928SBarry Smith 126ace3abfcSBarry 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) 1274b9ad928SBarry Smith { 128f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 129f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 130dfbe8321SBarry Smith PetscErrorCode ierr; 131391689abSStefano Zampini PC tpc; 132391689abSStefano Zampini PetscBool changeu,changed; 133f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 1344b9ad928SBarry Smith 1354b9ad928SBarry Smith PetscFunctionBegin; 136a762d673SBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 137a762d673SBarry Smith for (i=0; i<levels; i++) { 138a762d673SBarry Smith if (!mglevels[i]->A) { 139a762d673SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 140a762d673SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 141a762d673SBarry Smith } 142a762d673SBarry Smith } 143391689abSStefano Zampini 144391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 145391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 146391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 147391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 148391689abSStefano Zampini if (!changed && !changeu) { 149391689abSStefano Zampini ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 150f3fbd535SBarry Smith mglevels[levels-1]->b = b; 151391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 152391689abSStefano Zampini if (!mglevels[levels-1]->b) { 153391689abSStefano Zampini Vec *vec; 154391689abSStefano Zampini 155391689abSStefano Zampini ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 156391689abSStefano Zampini mglevels[levels-1]->b = *vec; 1577ae21d43SStefano Zampini ierr = PetscFree(vec);CHKERRQ(ierr); 158391689abSStefano Zampini } 159391689abSStefano Zampini ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 160391689abSStefano Zampini } 161f3fbd535SBarry Smith mglevels[levels-1]->x = x; 1624b9ad928SBarry Smith 16331567311SBarry Smith mg->rtol = rtol; 16431567311SBarry Smith mg->abstol = abstol; 16531567311SBarry Smith mg->dtol = dtol; 1664b9ad928SBarry Smith if (rtol) { 1674b9ad928SBarry Smith /* compute initial residual norm for relative convergence test */ 1684b9ad928SBarry Smith PetscReal rnorm; 1697319c654SBarry Smith if (zeroguess) { 1707319c654SBarry Smith ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); 1717319c654SBarry Smith } else { 172f3fbd535SBarry Smith ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); 1734b9ad928SBarry Smith ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); 1747319c654SBarry Smith } 17531567311SBarry Smith mg->ttol = PetscMax(rtol*rnorm,abstol); 1762fa5cd67SKarl Rupp } else if (abstol) mg->ttol = abstol; 1772fa5cd67SKarl Rupp else mg->ttol = 0.0; 1784b9ad928SBarry Smith 1794d0a8057SBarry Smith /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't 180336babb1SJed Brown stop prematurely due to small residual */ 1814d0a8057SBarry Smith for (i=1; i<levels; i++) { 182f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 183f3fbd535SBarry Smith if (mglevels[i]->smoothu != mglevels[i]->smoothd) { 18423067569SBarry Smith /* For Richardson the initial guess is nonzero since it is solving in each cycle the original system not just applying as a preconditioner */ 18523067569SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 186f3fbd535SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 1874b9ad928SBarry Smith } 1884d0a8057SBarry Smith } 1894d0a8057SBarry Smith 1904d0a8057SBarry Smith *reason = (PCRichardsonConvergedReason)0; 1914d0a8057SBarry Smith for (i=0; i<its; i++) { 19230b0564aSStefano Zampini ierr = PCMGMCycle_Private(pc,mglevels+levels-1,PETSC_FALSE,PETSC_FALSE,reason);CHKERRQ(ierr); 1934d0a8057SBarry Smith if (*reason) break; 1944d0a8057SBarry Smith } 1954d0a8057SBarry Smith if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; 1964d0a8057SBarry Smith *outits = i; 197391689abSStefano Zampini if (!changed && !changeu) mglevels[levels-1]->b = NULL; 1984b9ad928SBarry Smith PetscFunctionReturn(0); 1994b9ad928SBarry Smith } 2004b9ad928SBarry Smith 2013aeaf226SBarry Smith PetscErrorCode PCReset_MG(PC pc) 2023aeaf226SBarry Smith { 2033aeaf226SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 2043aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 2053aeaf226SBarry Smith PetscErrorCode ierr; 206f3b08a26SMatthew G. Knepley PetscInt i,c,n; 2073aeaf226SBarry Smith 2083aeaf226SBarry Smith PetscFunctionBegin; 2093aeaf226SBarry Smith if (mglevels) { 2103aeaf226SBarry Smith n = mglevels[0]->levels; 2113aeaf226SBarry Smith for (i=0; i<n-1; i++) { 2123aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i+1]->r);CHKERRQ(ierr); 2133aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->b);CHKERRQ(ierr); 2143aeaf226SBarry Smith ierr = VecDestroy(&mglevels[i]->x);CHKERRQ(ierr); 21530b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i+1]->R);CHKERRQ(ierr); 21630b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i]->B);CHKERRQ(ierr); 21730b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i]->X);CHKERRQ(ierr); 21841b6fd38SMatthew G. Knepley ierr = VecDestroy(&mglevels[i]->crx);CHKERRQ(ierr); 21941b6fd38SMatthew G. Knepley ierr = VecDestroy(&mglevels[i]->crb);CHKERRQ(ierr); 2203aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->restrct);CHKERRQ(ierr); 2213aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i+1]->interpolate);CHKERRQ(ierr); 2220de8493bSLawrence Mitchell ierr = MatDestroy(&mglevels[i+1]->inject);CHKERRQ(ierr); 22373250ac0SBarry Smith ierr = VecDestroy(&mglevels[i+1]->rscale);CHKERRQ(ierr); 2243aeaf226SBarry Smith } 22541b6fd38SMatthew G. Knepley ierr = VecDestroy(&mglevels[n-1]->crx);CHKERRQ(ierr); 22641b6fd38SMatthew G. Knepley ierr = VecDestroy(&mglevels[n-1]->crb);CHKERRQ(ierr); 227391689abSStefano Zampini /* this is not null only if the smoother on the finest level 228391689abSStefano Zampini changes the rhs during PreSolve */ 229391689abSStefano Zampini ierr = VecDestroy(&mglevels[n-1]->b);CHKERRQ(ierr); 23030b0564aSStefano Zampini ierr = MatDestroy(&mglevels[n-1]->B);CHKERRQ(ierr); 2313aeaf226SBarry Smith 2323aeaf226SBarry Smith for (i=0; i<n; i++) { 233f3b08a26SMatthew G. Knepley if (mglevels[i]->coarseSpace) for (c = 0; c < mg->Nc; ++c) {ierr = VecDestroy(&mglevels[i]->coarseSpace[c]);CHKERRQ(ierr);} 234f3b08a26SMatthew G. Knepley ierr = PetscFree(mglevels[i]->coarseSpace);CHKERRQ(ierr); 235f3b08a26SMatthew G. Knepley mglevels[i]->coarseSpace = NULL; 2363aeaf226SBarry Smith ierr = MatDestroy(&mglevels[i]->A);CHKERRQ(ierr); 2373aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 2383aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothd);CHKERRQ(ierr); 2393aeaf226SBarry Smith } 2403aeaf226SBarry Smith ierr = KSPReset(mglevels[i]->smoothu);CHKERRQ(ierr); 24141b6fd38SMatthew G. Knepley if (mglevels[i]->cr) {ierr = KSPReset(mglevels[i]->cr);CHKERRQ(ierr);} 2423aeaf226SBarry Smith } 243f3b08a26SMatthew G. Knepley mg->Nc = 0; 2443aeaf226SBarry Smith } 2453aeaf226SBarry Smith PetscFunctionReturn(0); 2463aeaf226SBarry Smith } 2473aeaf226SBarry Smith 24841b6fd38SMatthew G. Knepley /* Implementing CR 24941b6fd38SMatthew G. Knepley 25041b6fd38SMatthew 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 25141b6fd38SMatthew G. Knepley 25241b6fd38SMatthew G. Knepley Inj^T (Inj Inj^T)^{-1} Inj 25341b6fd38SMatthew G. Knepley 25441b6fd38SMatthew G. Knepley and if Inj is a VecScatter, as it is now in PETSc, we have 25541b6fd38SMatthew G. Knepley 25641b6fd38SMatthew G. Knepley Inj^T Inj 25741b6fd38SMatthew G. Knepley 25841b6fd38SMatthew G. Knepley and 25941b6fd38SMatthew G. Knepley 26041b6fd38SMatthew G. Knepley S = I - Inj^T Inj 26141b6fd38SMatthew G. Knepley 26241b6fd38SMatthew G. Knepley since 26341b6fd38SMatthew G. Knepley 26441b6fd38SMatthew G. Knepley Inj S = Inj - (Inj Inj^T) Inj = 0. 26541b6fd38SMatthew G. Knepley 26641b6fd38SMatthew G. Knepley Brannick suggests 26741b6fd38SMatthew G. Knepley 26841b6fd38SMatthew G. Knepley A \to S^T A S \qquad\mathrm{and}\qquad M \to S^T M S 26941b6fd38SMatthew G. Knepley 27041b6fd38SMatthew 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 27141b6fd38SMatthew G. Knepley 27241b6fd38SMatthew G. Knepley M^{-1} A \to S M^{-1} A S 27341b6fd38SMatthew G. Knepley 27441b6fd38SMatthew G. Knepley In fact, since it is somewhat hard in PETSc to do the symmetric application, we will just apply S on the left. 27541b6fd38SMatthew G. Knepley 27641b6fd38SMatthew G. Knepley Check: || Inj P - I ||_F < tol 27741b6fd38SMatthew G. Knepley Check: In general, Inj Inj^T = I 27841b6fd38SMatthew G. Knepley */ 27941b6fd38SMatthew G. Knepley 28041b6fd38SMatthew G. Knepley typedef struct { 28141b6fd38SMatthew G. Knepley PC mg; /* The PCMG object */ 28241b6fd38SMatthew G. Knepley PetscInt l; /* The multigrid level for this solver */ 28341b6fd38SMatthew G. Knepley Mat Inj; /* The injection matrix */ 28441b6fd38SMatthew G. Knepley Mat S; /* I - Inj^T Inj */ 28541b6fd38SMatthew G. Knepley } CRContext; 28641b6fd38SMatthew G. Knepley 28741b6fd38SMatthew G. Knepley static PetscErrorCode CRSetup_Private(PC pc) 28841b6fd38SMatthew G. Knepley { 28941b6fd38SMatthew G. Knepley CRContext *ctx; 29041b6fd38SMatthew G. Knepley Mat It; 29141b6fd38SMatthew G. Knepley PetscErrorCode ierr; 29241b6fd38SMatthew G. Knepley 29341b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 29441b6fd38SMatthew G. Knepley ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr); 29541b6fd38SMatthew G. Knepley ierr = PCMGGetInjection(ctx->mg, ctx->l, &It);CHKERRQ(ierr); 29641b6fd38SMatthew G. Knepley if (!It) SETERRQ(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONGSTATE, "CR requires that injection be defined for this PCMG"); 29741b6fd38SMatthew G. Knepley ierr = MatCreateTranspose(It, &ctx->Inj);CHKERRQ(ierr); 29841b6fd38SMatthew G. Knepley ierr = MatCreateNormal(ctx->Inj, &ctx->S);CHKERRQ(ierr); 29941b6fd38SMatthew G. Knepley ierr = MatScale(ctx->S, -1.0);CHKERRQ(ierr); 30041b6fd38SMatthew G. Knepley ierr = MatShift(ctx->S, 1.0);CHKERRQ(ierr); 30141b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 30241b6fd38SMatthew G. Knepley } 30341b6fd38SMatthew G. Knepley 30441b6fd38SMatthew G. Knepley static PetscErrorCode CRApply_Private(PC pc, Vec x, Vec y) 30541b6fd38SMatthew G. Knepley { 30641b6fd38SMatthew G. Knepley CRContext *ctx; 30741b6fd38SMatthew G. Knepley PetscErrorCode ierr; 30841b6fd38SMatthew G. Knepley 30941b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 31041b6fd38SMatthew G. Knepley ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr); 31141b6fd38SMatthew G. Knepley ierr = MatMult(ctx->S, x, y);CHKERRQ(ierr); 31241b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 31341b6fd38SMatthew G. Knepley } 31441b6fd38SMatthew G. Knepley 31541b6fd38SMatthew G. Knepley static PetscErrorCode CRDestroy_Private(PC pc) 31641b6fd38SMatthew G. Knepley { 31741b6fd38SMatthew G. Knepley CRContext *ctx; 31841b6fd38SMatthew G. Knepley PetscErrorCode ierr; 31941b6fd38SMatthew G. Knepley 32041b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 32141b6fd38SMatthew G. Knepley ierr = PCShellGetContext(pc, (void **) &ctx);CHKERRQ(ierr); 32241b6fd38SMatthew G. Knepley ierr = MatDestroy(&ctx->Inj);CHKERRQ(ierr); 32341b6fd38SMatthew G. Knepley ierr = MatDestroy(&ctx->S);CHKERRQ(ierr); 32441b6fd38SMatthew G. Knepley ierr = PetscFree(ctx);CHKERRQ(ierr); 32541b6fd38SMatthew G. Knepley ierr = PCShellSetContext(pc, NULL);CHKERRQ(ierr); 32641b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 32741b6fd38SMatthew G. Knepley } 32841b6fd38SMatthew G. Knepley 32941b6fd38SMatthew G. Knepley static PetscErrorCode CreateCR_Private(PC pc, PetscInt l, PC *cr) 33041b6fd38SMatthew G. Knepley { 33141b6fd38SMatthew G. Knepley CRContext *ctx; 33241b6fd38SMatthew G. Knepley PetscErrorCode ierr; 33341b6fd38SMatthew G. Knepley 33441b6fd38SMatthew G. Knepley PetscFunctionBeginUser; 33541b6fd38SMatthew G. Knepley ierr = PCCreate(PetscObjectComm((PetscObject) pc), cr);CHKERRQ(ierr); 33641b6fd38SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *cr, "S (complementary projector to injection)");CHKERRQ(ierr); 33741b6fd38SMatthew G. Knepley ierr = PetscCalloc1(1, &ctx);CHKERRQ(ierr); 33841b6fd38SMatthew G. Knepley ctx->mg = pc; 33941b6fd38SMatthew G. Knepley ctx->l = l; 34041b6fd38SMatthew G. Knepley ierr = PCSetType(*cr, PCSHELL);CHKERRQ(ierr); 34141b6fd38SMatthew G. Knepley ierr = PCShellSetContext(*cr, ctx);CHKERRQ(ierr); 34241b6fd38SMatthew G. Knepley ierr = PCShellSetApply(*cr, CRApply_Private);CHKERRQ(ierr); 34341b6fd38SMatthew G. Knepley ierr = PCShellSetSetUp(*cr, CRSetup_Private);CHKERRQ(ierr); 34441b6fd38SMatthew G. Knepley ierr = PCShellSetDestroy(*cr, CRDestroy_Private);CHKERRQ(ierr); 34541b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 34641b6fd38SMatthew G. Knepley } 34741b6fd38SMatthew G. Knepley 34897d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels_MG(PC pc,PetscInt levels,MPI_Comm *comms) 3494b9ad928SBarry Smith { 350dfbe8321SBarry Smith PetscErrorCode ierr; 351f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 352ce94432eSBarry Smith MPI_Comm comm; 3533aeaf226SBarry Smith PC_MG_Levels **mglevels = mg->levels; 35410eca3edSLisandro Dalcin PCMGType mgtype = mg->am; 35510167fecSBarry Smith PetscInt mgctype = (PetscInt) PC_MG_CYCLE_V; 356f3fbd535SBarry Smith PetscInt i; 357f3fbd535SBarry Smith PetscMPIInt size; 358f3fbd535SBarry Smith const char *prefix; 359f3fbd535SBarry Smith PC ipc; 3603aeaf226SBarry Smith PetscInt n; 3614b9ad928SBarry Smith 3624b9ad928SBarry Smith PetscFunctionBegin; 3630700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 364548767e0SJed Brown PetscValidLogicalCollectiveInt(pc,levels,2); 365548767e0SJed Brown if (mg->nlevels == levels) PetscFunctionReturn(0); 3661a97d7b6SStefano Zampini ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 3673aeaf226SBarry Smith if (mglevels) { 36810eca3edSLisandro Dalcin mgctype = mglevels[0]->cycles; 3693aeaf226SBarry Smith /* changing the number of levels so free up the previous stuff */ 3703aeaf226SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 3713aeaf226SBarry Smith n = mglevels[0]->levels; 3723aeaf226SBarry Smith for (i=0; i<n; i++) { 3733aeaf226SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 3743aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 3753aeaf226SBarry Smith } 3763aeaf226SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 37741b6fd38SMatthew G. Knepley ierr = KSPDestroy(&mglevels[i]->cr);CHKERRQ(ierr); 3783aeaf226SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 3793aeaf226SBarry Smith } 3803aeaf226SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 3813aeaf226SBarry Smith } 382f3fbd535SBarry Smith 383f3fbd535SBarry Smith mg->nlevels = levels; 384f3fbd535SBarry Smith 385785e854fSJed Brown ierr = PetscMalloc1(levels,&mglevels);CHKERRQ(ierr); 3863bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,levels*(sizeof(PC_MG*)));CHKERRQ(ierr); 387f3fbd535SBarry Smith 388f3fbd535SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 389f3fbd535SBarry Smith 390a9db87a2SMatthew G Knepley mg->stageApply = 0; 391f3fbd535SBarry Smith for (i=0; i<levels; i++) { 392b00a9115SJed Brown ierr = PetscNewLog(pc,&mglevels[i]);CHKERRQ(ierr); 3932fa5cd67SKarl Rupp 394f3fbd535SBarry Smith mglevels[i]->level = i; 395f3fbd535SBarry Smith mglevels[i]->levels = levels; 39610eca3edSLisandro Dalcin mglevels[i]->cycles = mgctype; 397336babb1SJed Brown mg->default_smoothu = 2; 398336babb1SJed Brown mg->default_smoothd = 2; 39963e6d426SJed Brown mglevels[i]->eventsmoothsetup = 0; 40063e6d426SJed Brown mglevels[i]->eventsmoothsolve = 0; 40163e6d426SJed Brown mglevels[i]->eventresidual = 0; 40263e6d426SJed Brown mglevels[i]->eventinterprestrict = 0; 403f3fbd535SBarry Smith 404f3fbd535SBarry Smith if (comms) comm = comms[i]; 405d5a8f7e6SBarry Smith if (comm != MPI_COMM_NULL) { 406f3fbd535SBarry Smith ierr = KSPCreate(comm,&mglevels[i]->smoothd);CHKERRQ(ierr); 407422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(mglevels[i]->smoothd,pc->erroriffailure);CHKERRQ(ierr); 4080ee9a668SBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)mglevels[i]->smoothd,(PetscObject)pc,levels-i);CHKERRQ(ierr); 4090ee9a668SBarry Smith ierr = KSPSetOptionsPrefix(mglevels[i]->smoothd,prefix);CHKERRQ(ierr); 4100ee9a668SBarry Smith ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->smoothd, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 4110ee9a668SBarry Smith if (i || levels == 1) { 4120ee9a668SBarry Smith char tprefix[128]; 4130ee9a668SBarry Smith 414336babb1SJed Brown ierr = KSPSetType(mglevels[i]->smoothd,KSPCHEBYSHEV);CHKERRQ(ierr); 4150059c7bdSJed Brown ierr = KSPSetConvergenceTest(mglevels[i]->smoothd,KSPConvergedSkip,NULL,NULL);CHKERRQ(ierr); 416336babb1SJed Brown ierr = KSPSetNormType(mglevels[i]->smoothd,KSP_NORM_NONE);CHKERRQ(ierr); 417336babb1SJed Brown ierr = KSPGetPC(mglevels[i]->smoothd,&ipc);CHKERRQ(ierr); 418336babb1SJed Brown ierr = PCSetType(ipc,PCSOR);CHKERRQ(ierr); 4190ee9a668SBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 420f3fbd535SBarry Smith 42141b6fd38SMatthew G. Knepley ierr = PetscSNPrintf(tprefix,128,"mg_levels_%d_",(int)i);CHKERRQ(ierr); 4220ee9a668SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[i]->smoothd,tprefix);CHKERRQ(ierr); 4230ee9a668SBarry Smith } else { 424f3fbd535SBarry Smith ierr = KSPAppendOptionsPrefix(mglevels[0]->smoothd,"mg_coarse_");CHKERRQ(ierr); 425f3fbd535SBarry Smith 4269dbfc187SHong Zhang /* coarse solve is (redundant) LU by default; set shifttype NONZERO to avoid annoying zero-pivot in LU preconditioner */ 427f3fbd535SBarry Smith ierr = KSPSetType(mglevels[0]->smoothd,KSPPREONLY);CHKERRQ(ierr); 428f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&ipc);CHKERRQ(ierr); 429ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 430f3fbd535SBarry Smith if (size > 1) { 431f3fbd535SBarry Smith ierr = PCSetType(ipc,PCREDUNDANT);CHKERRQ(ierr); 432f3fbd535SBarry Smith } else { 433f3fbd535SBarry Smith ierr = PCSetType(ipc,PCLU);CHKERRQ(ierr); 434f3fbd535SBarry Smith } 435753b7fb9SBarry Smith ierr = PCFactorSetShiftType(ipc,MAT_SHIFT_INBLOCKS);CHKERRQ(ierr); 436f3fbd535SBarry Smith } 4373bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)mglevels[i]->smoothd);CHKERRQ(ierr); 438d5a8f7e6SBarry Smith } 439f3fbd535SBarry Smith mglevels[i]->smoothu = mglevels[i]->smoothd; 44031567311SBarry Smith mg->rtol = 0.0; 44131567311SBarry Smith mg->abstol = 0.0; 44231567311SBarry Smith mg->dtol = 0.0; 44331567311SBarry Smith mg->ttol = 0.0; 44431567311SBarry Smith mg->cyclesperpcapply = 1; 445f3fbd535SBarry Smith } 446f3fbd535SBarry Smith mg->levels = mglevels; 44710eca3edSLisandro Dalcin ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 4484b9ad928SBarry Smith PetscFunctionReturn(0); 4494b9ad928SBarry Smith } 4504b9ad928SBarry Smith 45197d33e41SMatthew G. Knepley /*@C 45297d33e41SMatthew G. Knepley PCMGSetLevels - Sets the number of levels to use with MG. 45397d33e41SMatthew G. Knepley Must be called before any other MG routine. 45497d33e41SMatthew G. Knepley 45597d33e41SMatthew G. Knepley Logically Collective on PC 45697d33e41SMatthew G. Knepley 45797d33e41SMatthew G. Knepley Input Parameters: 45897d33e41SMatthew G. Knepley + pc - the preconditioner context 45997d33e41SMatthew G. Knepley . levels - the number of levels 46097d33e41SMatthew G. Knepley - comms - optional communicators for each level; this is to allow solving the coarser problems 461d5a8f7e6SBarry Smith on smaller sets of processes. For processes that are not included in the computation 462*05552d4bSJunchao Zhang you must pass MPI_COMM_NULL. Use comms = NULL to specify that all processes 463*05552d4bSJunchao Zhang should participate in each level of problem. 46497d33e41SMatthew G. Knepley 46597d33e41SMatthew G. Knepley Level: intermediate 46697d33e41SMatthew G. Knepley 46797d33e41SMatthew G. Knepley Notes: 46897d33e41SMatthew G. Knepley If the number of levels is one then the multigrid uses the -mg_levels prefix 46997d33e41SMatthew G. Knepley for setting the level options rather than the -mg_coarse prefix. 47097d33e41SMatthew G. Knepley 471d5a8f7e6SBarry Smith You can free the information in comms after this routine is called. 472d5a8f7e6SBarry Smith 473d5a8f7e6SBarry Smith The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level 474d5a8f7e6SBarry Smith are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on 475d5a8f7e6SBarry Smith the two levels, rank 0 in the original communicator will pass in an array of 2 communicators 476d5a8f7e6SBarry Smith of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators 477d5a8f7e6SBarry Smith the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate 478d5a8f7e6SBarry Smith in the coarse grid solve. 479d5a8f7e6SBarry Smith 480d5a8f7e6SBarry Smith Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one 481d5a8f7e6SBarry Smith must take special care in providing the restriction and interpolation operation. We recommend 482d5a8f7e6SBarry Smith providing these as two step operations; first perform a standard restriction or interpolation on 483d5a8f7e6SBarry Smith the full number of ranks for that level and then use an MPI call to copy the resulting vector 484*05552d4bSJunchao Zhang array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, note in both 485d5a8f7e6SBarry Smith cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and 486d5a8f7e6SBarry Smith recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries. 487d5a8f7e6SBarry Smith 488*05552d4bSJunchao Zhang Fortran Notes: 489*05552d4bSJunchao Zhang Use comms = PETSC_NULL_MPI_COMM as the equivalent of NULL in the C interface. Note PETSC_NULL_MPI_COMM 490*05552d4bSJunchao Zhang is not MPI_COMM_NULL. It is more like PETSC_NULL_INTEGER, PETSC_NULL_REAL etc. 491d5a8f7e6SBarry Smith 49297d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels() 49397d33e41SMatthew G. Knepley @*/ 49497d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 49597d33e41SMatthew G. Knepley { 49697d33e41SMatthew G. Knepley PetscErrorCode ierr; 49797d33e41SMatthew G. Knepley 49897d33e41SMatthew G. Knepley PetscFunctionBegin; 49997d33e41SMatthew G. Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 50097d33e41SMatthew G. Knepley if (comms) PetscValidPointer(comms,3); 50137b5128cSMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr); 50297d33e41SMatthew G. Knepley PetscFunctionReturn(0); 50397d33e41SMatthew G. Knepley } 50497d33e41SMatthew G. Knepley 505c07bf074SBarry Smith 506c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 507c07bf074SBarry Smith { 508c07bf074SBarry Smith PetscErrorCode ierr; 509a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 510a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 511a06653b4SBarry Smith PetscInt i,n; 512c07bf074SBarry Smith 513c07bf074SBarry Smith PetscFunctionBegin; 514a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 515a06653b4SBarry Smith if (mglevels) { 516a06653b4SBarry Smith n = mglevels[0]->levels; 517a06653b4SBarry Smith for (i=0; i<n; i++) { 518a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 5196bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 520a06653b4SBarry Smith } 5216bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 52241b6fd38SMatthew G. Knepley ierr = KSPDestroy(&mglevels[i]->cr);CHKERRQ(ierr); 523a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 524a06653b4SBarry Smith } 525a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 526a06653b4SBarry Smith } 527c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 528fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);CHKERRQ(ierr); 529fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);CHKERRQ(ierr); 530f3fbd535SBarry Smith PetscFunctionReturn(0); 531f3fbd535SBarry Smith } 532f3fbd535SBarry Smith 533f3fbd535SBarry Smith 534f3fbd535SBarry Smith /* 535f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 536f3fbd535SBarry Smith or full cycle of multigrid. 537f3fbd535SBarry Smith 538f3fbd535SBarry Smith Note: 539f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 540f3fbd535SBarry Smith */ 54130b0564aSStefano Zampini static PetscErrorCode PCApply_MG_Internal(PC pc,Vec b,Vec x,Mat B,Mat X,PetscBool transpose) 542f3fbd535SBarry Smith { 543f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 544f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 545f3fbd535SBarry Smith PetscErrorCode ierr; 546391689abSStefano Zampini PC tpc; 547f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 54830b0564aSStefano Zampini PetscBool changeu,changed,matapp; 549f3fbd535SBarry Smith 550f3fbd535SBarry Smith PetscFunctionBegin; 55130b0564aSStefano Zampini matapp = (PetscBool)(B && X); 552a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 553e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 554298cc208SBarry Smith for (i=0; i<levels; i++) { 555e1d8e5deSBarry Smith if (!mglevels[i]->A) { 55623ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 557298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 558e1d8e5deSBarry Smith } 559e1d8e5deSBarry Smith } 560e1d8e5deSBarry Smith 561391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 562391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 563391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 564391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 565391689abSStefano Zampini if (!changeu && !changed) { 56630b0564aSStefano Zampini if (matapp) { 56730b0564aSStefano Zampini ierr = MatDestroy(&mglevels[levels-1]->B);CHKERRQ(ierr); 56830b0564aSStefano Zampini mglevels[levels-1]->B = B; 56930b0564aSStefano Zampini } else { 570391689abSStefano Zampini ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 571f3fbd535SBarry Smith mglevels[levels-1]->b = b; 57230b0564aSStefano Zampini } 573391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 57430b0564aSStefano Zampini if (matapp) { 57530b0564aSStefano Zampini if (mglevels[levels-1]->B) { 57630b0564aSStefano Zampini PetscInt N1,N2; 57730b0564aSStefano Zampini PetscBool flg; 57830b0564aSStefano Zampini 57930b0564aSStefano Zampini ierr = MatGetSize(mglevels[levels-1]->B,NULL,&N1);CHKERRQ(ierr); 58030b0564aSStefano Zampini ierr = MatGetSize(B,NULL,&N2);CHKERRQ(ierr); 58130b0564aSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)mglevels[levels-1]->B,((PetscObject)B)->type_name,&flg);CHKERRQ(ierr); 58230b0564aSStefano Zampini if (N1 != N2 || !flg) { 58330b0564aSStefano Zampini ierr = MatDestroy(&mglevels[levels-1]->B);CHKERRQ(ierr); 58430b0564aSStefano Zampini } 58530b0564aSStefano Zampini } 58630b0564aSStefano Zampini if (!mglevels[levels-1]->B) { 58730b0564aSStefano Zampini ierr = MatDuplicate(B,MAT_COPY_VALUES,&mglevels[levels-1]->B);CHKERRQ(ierr); 58830b0564aSStefano Zampini } else { 58930b0564aSStefano Zampini ierr = MatCopy(B,mglevels[levels-1]->B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 59030b0564aSStefano Zampini } 59130b0564aSStefano Zampini } else { 592391689abSStefano Zampini if (!mglevels[levels-1]->b) { 593391689abSStefano Zampini Vec *vec; 594391689abSStefano Zampini 595391689abSStefano Zampini ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 596391689abSStefano Zampini mglevels[levels-1]->b = *vec; 5977ae21d43SStefano Zampini ierr = PetscFree(vec);CHKERRQ(ierr); 598391689abSStefano Zampini } 599391689abSStefano Zampini ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 600391689abSStefano Zampini } 60130b0564aSStefano Zampini } 60230b0564aSStefano Zampini if (matapp) { mglevels[levels-1]->X = X; } 60330b0564aSStefano Zampini else { mglevels[levels-1]->x = x; } 60430b0564aSStefano Zampini 60530b0564aSStefano Zampini /* If coarser Xs are present, it means we have already block applied the PC at least once 60630b0564aSStefano Zampini Reset operators if sizes/type do no match */ 60730b0564aSStefano Zampini if (matapp && levels > 1 && mglevels[levels-2]->X) { 60830b0564aSStefano Zampini PetscInt Xc,Bc; 60930b0564aSStefano Zampini PetscBool flg; 61030b0564aSStefano Zampini 61130b0564aSStefano Zampini ierr = MatGetSize(mglevels[levels-2]->X,NULL,&Xc);CHKERRQ(ierr); 61230b0564aSStefano Zampini ierr = MatGetSize(mglevels[levels-1]->B,NULL,&Bc);CHKERRQ(ierr); 61330b0564aSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)mglevels[levels-2]->X,((PetscObject)mglevels[levels-1]->X)->type_name,&flg);CHKERRQ(ierr); 61430b0564aSStefano Zampini if (Xc != Bc || !flg) { 61530b0564aSStefano Zampini ierr = MatDestroy(&mglevels[levels-1]->R);CHKERRQ(ierr); 61630b0564aSStefano Zampini for (i=0;i<levels-1;i++) { 61730b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i]->R);CHKERRQ(ierr); 61830b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i]->B);CHKERRQ(ierr); 61930b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i]->X);CHKERRQ(ierr); 62030b0564aSStefano Zampini } 62130b0564aSStefano Zampini } 62230b0564aSStefano Zampini } 623391689abSStefano Zampini 62431567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 62530b0564aSStefano Zampini if (matapp) { ierr = MatZeroEntries(X);CHKERRQ(ierr); } 62630b0564aSStefano Zampini else { ierr = VecZeroEntries(x);CHKERRQ(ierr); } 62731567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 62830b0564aSStefano Zampini ierr = PCMGMCycle_Private(pc,mglevels+levels-1,transpose,matapp,NULL);CHKERRQ(ierr); 629f3fbd535SBarry Smith } 6302fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 63130b0564aSStefano Zampini ierr = PCMGACycle_Private(pc,mglevels,transpose,matapp);CHKERRQ(ierr); 6322fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 63330b0564aSStefano Zampini ierr = PCMGKCycle_Private(pc,mglevels,transpose,matapp);CHKERRQ(ierr); 6342fa5cd67SKarl Rupp } else { 63530b0564aSStefano Zampini ierr = PCMGFCycle_Private(pc,mglevels,transpose,matapp);CHKERRQ(ierr); 636f3fbd535SBarry Smith } 637a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 63830b0564aSStefano Zampini if (!changeu && !changed) { 63930b0564aSStefano Zampini if (matapp) { mglevels[levels-1]->B = NULL; } 64030b0564aSStefano Zampini else { mglevels[levels-1]->b = NULL; } 64130b0564aSStefano Zampini } 642f3fbd535SBarry Smith PetscFunctionReturn(0); 643f3fbd535SBarry Smith } 644f3fbd535SBarry Smith 645fcb023d4SJed Brown static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 646fcb023d4SJed Brown { 647fcb023d4SJed Brown PetscErrorCode ierr; 648fcb023d4SJed Brown 649fcb023d4SJed Brown PetscFunctionBegin; 65030b0564aSStefano Zampini ierr = PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_FALSE);CHKERRQ(ierr); 651fcb023d4SJed Brown PetscFunctionReturn(0); 652fcb023d4SJed Brown } 653fcb023d4SJed Brown 654fcb023d4SJed Brown static PetscErrorCode PCApplyTranspose_MG(PC pc,Vec b,Vec x) 655fcb023d4SJed Brown { 656fcb023d4SJed Brown PetscErrorCode ierr; 657fcb023d4SJed Brown 658fcb023d4SJed Brown PetscFunctionBegin; 65930b0564aSStefano Zampini ierr = PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_TRUE);CHKERRQ(ierr); 66030b0564aSStefano Zampini PetscFunctionReturn(0); 66130b0564aSStefano Zampini } 66230b0564aSStefano Zampini 66330b0564aSStefano Zampini static PetscErrorCode PCMatApply_MG(PC pc,Mat b,Mat x) 66430b0564aSStefano Zampini { 66530b0564aSStefano Zampini PetscErrorCode ierr; 66630b0564aSStefano Zampini 66730b0564aSStefano Zampini PetscFunctionBegin; 66830b0564aSStefano Zampini ierr = PCApply_MG_Internal(pc,NULL,NULL,b,x,PETSC_FALSE);CHKERRQ(ierr); 669fcb023d4SJed Brown PetscFunctionReturn(0); 670fcb023d4SJed Brown } 671f3fbd535SBarry Smith 6724416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 673f3fbd535SBarry Smith { 674f3fbd535SBarry Smith PetscErrorCode ierr; 675710315b6SLawrence Mitchell PetscInt levels,cycles; 676f3b08a26SMatthew G. Knepley PetscBool flg, flg2; 677f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 6783d3eaba7SBarry Smith PC_MG_Levels **mglevels; 679f3fbd535SBarry Smith PCMGType mgtype; 680f3fbd535SBarry Smith PCMGCycleType mgctype; 6812134b1e4SBarry Smith PCMGGalerkinType gtype; 682f3fbd535SBarry Smith 683f3fbd535SBarry Smith PetscFunctionBegin; 6841d743356SStefano Zampini levels = PetscMax(mg->nlevels,1); 685e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr); 686f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 6871a97d7b6SStefano Zampini if (!flg && !mg->levels && pc->dm) { 688eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 689eb3f98d2SBarry Smith levels++; 6903aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 691eb3f98d2SBarry Smith } 6920298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 6933aeaf226SBarry Smith mglevels = mg->levels; 6943aeaf226SBarry Smith 695f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 696f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 697f3fbd535SBarry Smith if (flg) { 698f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 6992fa5cd67SKarl Rupp } 7002134b1e4SBarry Smith gtype = mg->galerkin; 7012134b1e4SBarry Smith ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg);CHKERRQ(ierr); 7022134b1e4SBarry Smith if (flg) { 7032134b1e4SBarry Smith ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr); 704f3fbd535SBarry Smith } 705f3b08a26SMatthew G. Knepley flg2 = PETSC_FALSE; 706f3b08a26SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr); 707f3b08a26SMatthew G. Knepley if (flg) {ierr = PCMGSetAdaptInterpolation(pc, flg2);CHKERRQ(ierr);} 708f3b08a26SMatthew 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); 709f3b08a26SMatthew 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); 710f3b08a26SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);CHKERRQ(ierr); 71141b6fd38SMatthew G. Knepley flg2 = PETSC_FALSE; 71241b6fd38SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr); 71341b6fd38SMatthew G. Knepley if (flg) {ierr = PCMGSetAdaptCR(pc, flg2);CHKERRQ(ierr);} 714f56b1265SBarry Smith flg = PETSC_FALSE; 7158e5aa403SBarry Smith ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr); 716f442ab6aSBarry Smith if (flg) { 717f442ab6aSBarry Smith ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr); 718f442ab6aSBarry Smith } 71931567311SBarry Smith mgtype = mg->am; 720f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 721f3fbd535SBarry Smith if (flg) { 722f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 723f3fbd535SBarry Smith } 72431567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 7255363c1e0SLisandro Dalcin ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 726f3fbd535SBarry Smith if (flg) { 727f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 728f3fbd535SBarry Smith } 729f3fbd535SBarry Smith } 730f3fbd535SBarry Smith flg = PETSC_FALSE; 7310298fd71SBarry Smith ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 732f3fbd535SBarry Smith if (flg) { 733f3fbd535SBarry Smith PetscInt i; 734f3fbd535SBarry Smith char eventname[128]; 7351a97d7b6SStefano Zampini 736f3fbd535SBarry Smith levels = mglevels[0]->levels; 737f3fbd535SBarry Smith for (i=0; i<levels; i++) { 738f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 73963e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 740f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 74163e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 742f3fbd535SBarry Smith if (i) { 743f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 74463e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 745f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 74663e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 747f3fbd535SBarry Smith } 748f3fbd535SBarry Smith } 749eec35531SMatthew G Knepley 750ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 751eec35531SMatthew G Knepley { 752eec35531SMatthew G Knepley const char *sname = "MG Apply"; 753eec35531SMatthew G Knepley PetscStageLog stageLog; 754eec35531SMatthew G Knepley PetscInt st; 755eec35531SMatthew G Knepley 756eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 757eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 758eec35531SMatthew G Knepley PetscBool same; 759eec35531SMatthew G Knepley 760eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 7612fa5cd67SKarl Rupp if (same) mg->stageApply = st; 762eec35531SMatthew G Knepley } 763eec35531SMatthew G Knepley if (!mg->stageApply) { 764eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 765eec35531SMatthew G Knepley } 766eec35531SMatthew G Knepley } 767ec60ca73SSatish Balay #endif 768f3fbd535SBarry Smith } 769f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 770f3b08a26SMatthew G. Knepley /* Check option consistency */ 771f3b08a26SMatthew G. Knepley ierr = PCMGGetGalerkin(pc, >ype);CHKERRQ(ierr); 772f3b08a26SMatthew G. Knepley ierr = PCMGGetAdaptInterpolation(pc, &flg);CHKERRQ(ierr); 773f3b08a26SMatthew 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"); 774f3fbd535SBarry Smith PetscFunctionReturn(0); 775f3fbd535SBarry Smith } 776f3fbd535SBarry Smith 7770a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL}; 7780a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL}; 7790a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL}; 780f3b08a26SMatthew G. Knepley const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL}; 781f3fbd535SBarry Smith 7829804daf3SBarry Smith #include <petscdraw.h> 783f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 784f3fbd535SBarry Smith { 785f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 786f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 787f3fbd535SBarry Smith PetscErrorCode ierr; 788e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 789219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 790f3fbd535SBarry Smith 791f3fbd535SBarry Smith PetscFunctionBegin; 792251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 7935b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 794219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 795f3fbd535SBarry Smith if (iascii) { 796e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 797efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 79831567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 79931567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 800f3fbd535SBarry Smith } 8012134b1e4SBarry Smith if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 802f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 8032134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 8042134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr); 8052134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 8062134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr); 8072134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 8082134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr); 8094f66f45eSBarry Smith } else { 8104f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 811f3fbd535SBarry Smith } 8125adeb434SBarry Smith if (mg->view){ 8135adeb434SBarry Smith ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr); 8145adeb434SBarry Smith } 815f3fbd535SBarry Smith for (i=0; i<levels; i++) { 816f3fbd535SBarry Smith if (!i) { 81789cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 818f3fbd535SBarry Smith } else { 81989cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 820f3fbd535SBarry Smith } 821f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 822f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 823f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 824f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 825f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 826f3fbd535SBarry Smith } else if (i) { 82789cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 828f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 829f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 830f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 831f3fbd535SBarry Smith } 83241b6fd38SMatthew G. Knepley if (i && mglevels[i]->cr) { 83341b6fd38SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer,"CR solver on level %D -------------------------------\n",i);CHKERRQ(ierr); 83441b6fd38SMatthew G. Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 83541b6fd38SMatthew G. Knepley ierr = KSPView(mglevels[i]->cr,viewer);CHKERRQ(ierr); 83641b6fd38SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 83741b6fd38SMatthew G. Knepley } 838f3fbd535SBarry Smith } 8395b0b0462SBarry Smith } else if (isbinary) { 8405b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 8415b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 8425b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 8435b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 8445b0b0462SBarry Smith } 8455b0b0462SBarry Smith } 846219b1045SBarry Smith } else if (isdraw) { 847219b1045SBarry Smith PetscDraw draw; 848b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 849219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 850219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 8510298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 852219b1045SBarry Smith bottom = y - th; 853219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 854b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 855219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 856219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 857219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 858b4375e8dSPeter Brune } else { 859b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 860b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 861b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 862b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 863b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 864b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 865b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 866b4375e8dSPeter Brune } 8670298fd71SBarry Smith ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 8681cd381d6SBarry Smith bottom -= th; 869219b1045SBarry Smith } 8705b0b0462SBarry Smith } 871f3fbd535SBarry Smith PetscFunctionReturn(0); 872f3fbd535SBarry Smith } 873f3fbd535SBarry Smith 874af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 875cab2e9ccSBarry Smith 876f3fbd535SBarry Smith /* 877f3fbd535SBarry Smith Calls setup for the KSP on each level 878f3fbd535SBarry Smith */ 879f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 880f3fbd535SBarry Smith { 881f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 882f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 883f3fbd535SBarry Smith PetscErrorCode ierr; 8841a97d7b6SStefano Zampini PetscInt i,n; 88598e04984SBarry Smith PC cpc; 8863db492dfSBarry Smith PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 887f3fbd535SBarry Smith Mat dA,dB; 888f3fbd535SBarry Smith Vec tvec; 889218a07d4SBarry Smith DM *dms; 8900a545947SLisandro Dalcin PetscViewer viewer = NULL; 89141b6fd38SMatthew G. Knepley PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE; 892f3fbd535SBarry Smith 893f3fbd535SBarry Smith PetscFunctionBegin; 8941a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up"); 8951a97d7b6SStefano Zampini n = mglevels[0]->levels; 89601bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 8973aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 8983aeaf226SBarry Smith PetscInt levels; 8993aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 9003aeaf226SBarry Smith levels++; 9013aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 9020298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 9033aeaf226SBarry Smith n = levels; 9043aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 9053aeaf226SBarry Smith mglevels = mg->levels; 9063aeaf226SBarry Smith } 9073aeaf226SBarry Smith } 90898e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 9093aeaf226SBarry Smith 910f3fbd535SBarry Smith 911f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 912f3fbd535SBarry Smith /* so use those from global PC */ 913f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 9140298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 91598e04984SBarry Smith if (opsset) { 91698e04984SBarry Smith Mat mmat; 91723ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 91898e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 91998e04984SBarry Smith } 9204a5f07a5SMark F. Adams 92141b6fd38SMatthew G. Knepley /* Create CR solvers */ 92241b6fd38SMatthew G. Knepley ierr = PCMGGetAdaptCR(pc, &doCR);CHKERRQ(ierr); 92341b6fd38SMatthew G. Knepley if (doCR) { 92441b6fd38SMatthew G. Knepley const char *prefix; 92541b6fd38SMatthew G. Knepley 92641b6fd38SMatthew G. Knepley ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 92741b6fd38SMatthew G. Knepley for (i = 1; i < n; ++i) { 92841b6fd38SMatthew G. Knepley PC ipc, cr; 92941b6fd38SMatthew G. Knepley char crprefix[128]; 93041b6fd38SMatthew G. Knepley 93141b6fd38SMatthew G. Knepley ierr = KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr);CHKERRQ(ierr); 93241b6fd38SMatthew G. Knepley ierr = KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE);CHKERRQ(ierr); 93341b6fd38SMatthew G. Knepley ierr = PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i);CHKERRQ(ierr); 93441b6fd38SMatthew G. Knepley ierr = KSPSetOptionsPrefix(mglevels[i]->cr, prefix);CHKERRQ(ierr); 93541b6fd38SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 93641b6fd38SMatthew G. Knepley ierr = KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV);CHKERRQ(ierr); 93741b6fd38SMatthew G. Knepley ierr = KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL);CHKERRQ(ierr); 93841b6fd38SMatthew G. Knepley ierr = KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED);CHKERRQ(ierr); 93941b6fd38SMatthew G. Knepley ierr = KSPGetPC(mglevels[i]->cr, &ipc);CHKERRQ(ierr); 94041b6fd38SMatthew G. Knepley 94141b6fd38SMatthew G. Knepley ierr = PCSetType(ipc, PCCOMPOSITE);CHKERRQ(ierr); 94241b6fd38SMatthew G. Knepley ierr = PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 94341b6fd38SMatthew G. Knepley ierr = PCCompositeAddPCType(ipc, PCSOR);CHKERRQ(ierr); 94441b6fd38SMatthew G. Knepley ierr = CreateCR_Private(pc, i, &cr);CHKERRQ(ierr); 94541b6fd38SMatthew G. Knepley ierr = PCCompositeAddPC(ipc, cr);CHKERRQ(ierr); 94641b6fd38SMatthew G. Knepley ierr = PCDestroy(&cr);CHKERRQ(ierr); 94741b6fd38SMatthew G. Knepley 94841b6fd38SMatthew G. Knepley ierr = KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 94941b6fd38SMatthew G. Knepley ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE);CHKERRQ(ierr); 95041b6fd38SMatthew G. Knepley ierr = PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i);CHKERRQ(ierr); 95141b6fd38SMatthew G. Knepley ierr = KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix);CHKERRQ(ierr); 95241b6fd38SMatthew G. Knepley ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr);CHKERRQ(ierr); 95341b6fd38SMatthew G. Knepley } 95441b6fd38SMatthew G. Knepley } 95541b6fd38SMatthew G. Knepley 9564a5f07a5SMark F. Adams if (!opsset) { 95771b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 9584a5f07a5SMark F. Adams if (use_amat) { 959f3fbd535SBarry 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); 96023ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 96169858f1bSStefano Zampini } else { 9624a5f07a5SMark 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); 96323ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 9644a5f07a5SMark F. Adams } 9654a5f07a5SMark F. Adams } 966f3fbd535SBarry Smith 9679c7ffb39SBarry Smith for (i=n-1; i>0; i--) { 9689c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 9699c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 9709c7ffb39SBarry Smith continue; 9719c7ffb39SBarry Smith } 9729c7ffb39SBarry Smith } 9732134b1e4SBarry Smith 9742134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 9752134b1e4SBarry Smith if (dA == dB) dAeqdB = PETSC_TRUE; 9762134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) { 9772134b1e4SBarry Smith needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 9782134b1e4SBarry Smith } 9792134b1e4SBarry Smith 9802134b1e4SBarry Smith 9819c7ffb39SBarry Smith /* 9829c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 9839c7ffb39SBarry Smith Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 9849c7ffb39SBarry Smith */ 9852134b1e4SBarry Smith if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 9862d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 9872e499ae9SBarry Smith Mat p; 98873250ac0SBarry Smith Vec rscale; 989785e854fSJed Brown ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 990218a07d4SBarry Smith dms[n-1] = pc->dm; 991ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 992ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 99341fce8fdSHong Zhang /* 99441fce8fdSHong Zhang Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level. 99541fce8fdSHong Zhang Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options. 99641fce8fdSHong Zhang But it is safe to use -dm_mat_type. 99741fce8fdSHong Zhang 99841fce8fdSHong Zhang The mat type should not be hardcoded like this, we need to find a better way. 99941fce8fdSHong Zhang ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr); 100041fce8fdSHong Zhang */ 1001218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 1002942e3340SBarry Smith DMKSP kdm; 1003eab52d2bSLawrence Mitchell PetscBool dmhasrestrict, dmhasinject; 10043c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 10052134b1e4SBarry Smith if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 1006c27ee7a3SPatrick Farrell if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1007c27ee7a3SPatrick Farrell ierr = KSPSetDM(mglevels[i]->smoothu,dms[i]);CHKERRQ(ierr); 1008c27ee7a3SPatrick Farrell if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);CHKERRQ(ierr);} 1009c27ee7a3SPatrick Farrell } 101041b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 101141b6fd38SMatthew G. Knepley ierr = KSPSetDM(mglevels[i]->cr,dms[i]);CHKERRQ(ierr); 101241b6fd38SMatthew G. Knepley if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE);CHKERRQ(ierr);} 101341b6fd38SMatthew G. Knepley } 1014942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 1015d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 1016d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 10170298fd71SBarry Smith kdm->ops->computerhs = NULL; 10180298fd71SBarry Smith kdm->rhsctx = NULL; 101924c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 1020e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 10212d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 102200ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 102373250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 10246bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 1025218a07d4SBarry Smith } 10263ad4599aSBarry Smith ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr); 10273ad4599aSBarry Smith if (dmhasrestrict && !mglevels[i+1]->restrct){ 10283ad4599aSBarry Smith ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr); 10293ad4599aSBarry Smith ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr); 10303ad4599aSBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 10313ad4599aSBarry Smith } 1032eab52d2bSLawrence Mitchell ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr); 1033eab52d2bSLawrence Mitchell if (dmhasinject && !mglevels[i+1]->inject){ 1034eab52d2bSLawrence Mitchell ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr); 1035eab52d2bSLawrence Mitchell ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr); 1036eab52d2bSLawrence Mitchell ierr = MatDestroy(&p);CHKERRQ(ierr); 1037eab52d2bSLawrence Mitchell } 103824c3aa18SBarry Smith } 10392d2b81a6SBarry Smith 1040ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 10412d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 1042ce4cda84SJed Brown } 1043cab2e9ccSBarry Smith 1044ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 10452134b1e4SBarry Smith /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 1046cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 1047cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 1048c27ee7a3SPatrick Farrell if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) { 1049c27ee7a3SPatrick Farrell ierr = KSPSetDM(mglevels[n-1]->smoothu,pc->dm);CHKERRQ(ierr); 1050c27ee7a3SPatrick Farrell ierr = KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);CHKERRQ(ierr); 1051c27ee7a3SPatrick Farrell } 105241b6fd38SMatthew G. Knepley if (mglevels[n-1]->cr) { 105341b6fd38SMatthew G. Knepley ierr = KSPSetDM(mglevels[n-1]->cr,pc->dm);CHKERRQ(ierr); 105441b6fd38SMatthew G. Knepley ierr = KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE);CHKERRQ(ierr); 105541b6fd38SMatthew G. Knepley } 1056218a07d4SBarry Smith } 1057218a07d4SBarry Smith 10582134b1e4SBarry Smith if (mg->galerkin < PC_MG_GALERKIN_NONE) { 10592134b1e4SBarry Smith Mat A,B; 10602134b1e4SBarry Smith PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 10612134b1e4SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 10622134b1e4SBarry Smith 10632134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 10642134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE; 10652134b1e4SBarry Smith if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 1066f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 1067b9367914SBarry 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"); 1068b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 1069b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 1070b9367914SBarry Smith } 1071b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 1072b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 1073b9367914SBarry Smith } 10742134b1e4SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 10752134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr); 10762134b1e4SBarry Smith } 10772134b1e4SBarry Smith if (doA) { 10782df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr); 10792134b1e4SBarry Smith } 10802134b1e4SBarry Smith if (doB) { 10812df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr); 1082b9367914SBarry Smith } 10832134b1e4SBarry Smith /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 10842134b1e4SBarry Smith if (!doA && dAeqdB) { 10852134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);} 10862134b1e4SBarry Smith A = B; 10872134b1e4SBarry Smith } else if (!doA && reuse == MAT_INITIAL_MATRIX) { 10882134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr); 10892134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1090b9367914SBarry Smith } 10912134b1e4SBarry Smith if (!doB && dAeqdB) { 10922134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);} 10932134b1e4SBarry Smith B = A; 10942134b1e4SBarry Smith } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 109523ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 10962134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 10977d537d92SJed Brown } 10982134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 10992134b1e4SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr); 11002134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr); 11012134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr); 11022134b1e4SBarry Smith } 11032134b1e4SBarry Smith dA = A; 1104f3fbd535SBarry Smith dB = B; 1105f3fbd535SBarry Smith } 1106f3fbd535SBarry Smith } 1107f3b08a26SMatthew G. Knepley 1108f3b08a26SMatthew G. Knepley 1109f3b08a26SMatthew G. Knepley /* Adapt interpolation matrices */ 1110f3b08a26SMatthew G. Knepley if (mg->adaptInterpolation) { 1111f3b08a26SMatthew G. Knepley mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */ 1112f3b08a26SMatthew G. Knepley for (i = 0; i < n; ++i) { 1113f3b08a26SMatthew G. Knepley ierr = PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);CHKERRQ(ierr); 1114f3b08a26SMatthew 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);} 1115f3b08a26SMatthew G. Knepley } 1116f3b08a26SMatthew G. Knepley for (i = n-2; i > -1; --i) { 1117f3b08a26SMatthew G. Knepley ierr = PCMGRecomputeLevelOperators_Internal(pc, i);CHKERRQ(ierr); 1118f3b08a26SMatthew G. Knepley } 1119f3b08a26SMatthew G. Knepley } 1120f3b08a26SMatthew G. Knepley 11212134b1e4SBarry Smith if (needRestricts && pc->dm) { 1122caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 1123ccceb50cSJed Brown DM dmfine,dmcoarse; 1124caa4e7f2SJed Brown Mat Restrict,Inject; 1125caa4e7f2SJed Brown Vec rscale; 1126ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 1127ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 1128caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 1129caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 1130eab52d2bSLawrence Mitchell ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr); 1131caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 1132caa4e7f2SJed Brown } 1133caa4e7f2SJed Brown } 1134f3fbd535SBarry Smith 1135f3fbd535SBarry Smith if (!pc->setupcalled) { 1136f3fbd535SBarry Smith for (i=0; i<n; i++) { 1137f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 1138f3fbd535SBarry Smith } 1139f3fbd535SBarry Smith for (i=1; i<n; i++) { 1140f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 1141f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 1142f3fbd535SBarry Smith } 114341b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 114441b6fd38SMatthew G. Knepley ierr = KSPSetFromOptions(mglevels[i]->cr);CHKERRQ(ierr); 114541b6fd38SMatthew G. Knepley } 1146f3fbd535SBarry Smith } 11473ad4599aSBarry Smith /* insure that if either interpolation or restriction is set the other other one is set */ 1148f3fbd535SBarry Smith for (i=1; i<n; i++) { 11493ad4599aSBarry Smith ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr); 11503ad4599aSBarry Smith ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr); 1151f3fbd535SBarry Smith } 1152f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 1153f3fbd535SBarry Smith if (!mglevels[i]->b) { 1154f3fbd535SBarry Smith Vec *vec; 11552a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 1156f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 11576bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 1158f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 1159f3fbd535SBarry Smith } 1160f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 1161f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 1162f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 11636bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 1164f3fbd535SBarry Smith } 1165f3fbd535SBarry Smith if (!mglevels[i]->x) { 1166f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 1167f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 11686bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 1169f3fbd535SBarry Smith } 117041b6fd38SMatthew G. Knepley if (doCR) { 117141b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crx);CHKERRQ(ierr); 117241b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crb);CHKERRQ(ierr); 117341b6fd38SMatthew G. Knepley ierr = VecZeroEntries(mglevels[i]->crb);CHKERRQ(ierr); 117441b6fd38SMatthew G. Knepley } 1175f3fbd535SBarry Smith } 1176f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 1177f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 1178f3fbd535SBarry Smith Vec *vec; 11792a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 1180f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 11816bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 1182f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 1183f3fbd535SBarry Smith } 118441b6fd38SMatthew G. Knepley if (doCR) { 118541b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx);CHKERRQ(ierr); 118641b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb);CHKERRQ(ierr); 118741b6fd38SMatthew G. Knepley ierr = VecZeroEntries(mglevels[n-1]->crb);CHKERRQ(ierr); 118841b6fd38SMatthew G. Knepley } 1189f3fbd535SBarry Smith } 1190f3fbd535SBarry Smith 119198e04984SBarry Smith if (pc->dm) { 119298e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 119398e04984SBarry Smith for (i=0; i<n-1; i++) { 119498e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 119598e04984SBarry Smith } 119698e04984SBarry Smith } 1197f3fbd535SBarry Smith 1198f3fbd535SBarry Smith for (i=1; i<n; i++) { 11992a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 1200f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 1201f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 1202f3fbd535SBarry Smith } 120341b6fd38SMatthew G. Knepley if (mglevels[i]->cr) {ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr);} 120463e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1205f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 1206c0decd05SBarry Smith if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1207899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1208899639b0SHong Zhang } 120963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1210d42688cbSBarry Smith if (!mglevels[i]->residual) { 1211d42688cbSBarry Smith Mat mat; 121213842ffbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 121354b2cd4bSJed Brown ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 1214d42688cbSBarry Smith } 1215fcb023d4SJed Brown if (!mglevels[i]->residualtranspose) { 1216fcb023d4SJed Brown Mat mat; 1217fcb023d4SJed Brown ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 1218fcb023d4SJed Brown ierr = PCMGSetResidualTranspose(pc,i,PCMGResidualTransposeDefault,mat);CHKERRQ(ierr); 1219fcb023d4SJed Brown } 1220f3fbd535SBarry Smith } 1221f3fbd535SBarry Smith for (i=1; i<n; i++) { 1222f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 1223f3fbd535SBarry Smith Mat downmat,downpmat; 1224f3fbd535SBarry Smith 1225f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 12260298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 1227f3fbd535SBarry Smith if (!opsset) { 122823ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 122923ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 1230f3fbd535SBarry Smith } 1231f3fbd535SBarry Smith 1232f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 123363e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1234f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 1235c0decd05SBarry Smith if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) { 1236899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1237899639b0SHong Zhang } 123863e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1239f3fbd535SBarry Smith } 124041b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 124141b6fd38SMatthew G. Knepley Mat downmat,downpmat; 124241b6fd38SMatthew G. Knepley 124341b6fd38SMatthew G. Knepley /* check if operators have been set for up, if not use down operators to set them */ 124441b6fd38SMatthew G. Knepley ierr = KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL);CHKERRQ(ierr); 124541b6fd38SMatthew G. Knepley if (!opsset) { 124641b6fd38SMatthew G. Knepley ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 124741b6fd38SMatthew G. Knepley ierr = KSPSetOperators(mglevels[i]->cr,downmat,downpmat);CHKERRQ(ierr); 124841b6fd38SMatthew G. Knepley } 124941b6fd38SMatthew G. Knepley 125041b6fd38SMatthew G. Knepley ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr); 125141b6fd38SMatthew G. Knepley if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 125241b6fd38SMatthew G. Knepley ierr = KSPSetUp(mglevels[i]->cr);CHKERRQ(ierr); 125341b6fd38SMatthew G. Knepley if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) { 125441b6fd38SMatthew G. Knepley pc->failedreason = PC_SUBPC_ERROR; 125541b6fd38SMatthew G. Knepley } 125641b6fd38SMatthew G. Knepley if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 125741b6fd38SMatthew G. Knepley } 1258f3fbd535SBarry Smith } 1259f3fbd535SBarry Smith 126063e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1261f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 1262c0decd05SBarry Smith if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1263899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1264899639b0SHong Zhang } 126563e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1266f3fbd535SBarry Smith 1267f3fbd535SBarry Smith /* 1268f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 1269e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 1270f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 1271f3fbd535SBarry Smith 1272f3fbd535SBarry Smith Only support one or the other at the same time. 1273f3fbd535SBarry Smith */ 1274f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 1275c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 1276ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 1277f3fbd535SBarry Smith dump = PETSC_FALSE; 1278f3fbd535SBarry Smith #endif 1279c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 1280ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 1281f3fbd535SBarry Smith 1282f3fbd535SBarry Smith if (viewer) { 1283f3fbd535SBarry Smith for (i=1; i<n; i++) { 1284f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 1285f3fbd535SBarry Smith } 1286f3fbd535SBarry Smith for (i=0; i<n; i++) { 1287f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 1288f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 1289f3fbd535SBarry Smith } 1290f3fbd535SBarry Smith } 1291f3fbd535SBarry Smith PetscFunctionReturn(0); 1292f3fbd535SBarry Smith } 1293f3fbd535SBarry Smith 1294f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 1295f3fbd535SBarry Smith 129697d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 129797d33e41SMatthew G. Knepley { 129897d33e41SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 129997d33e41SMatthew G. Knepley 130097d33e41SMatthew G. Knepley PetscFunctionBegin; 130197d33e41SMatthew G. Knepley *levels = mg->nlevels; 130297d33e41SMatthew G. Knepley PetscFunctionReturn(0); 130397d33e41SMatthew G. Knepley } 130497d33e41SMatthew G. Knepley 13054b9ad928SBarry Smith /*@ 130697177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 13074b9ad928SBarry Smith 13084b9ad928SBarry Smith Not Collective 13094b9ad928SBarry Smith 13104b9ad928SBarry Smith Input Parameter: 13114b9ad928SBarry Smith . pc - the preconditioner context 13124b9ad928SBarry Smith 13134b9ad928SBarry Smith Output parameter: 13144b9ad928SBarry Smith . levels - the number of levels 13154b9ad928SBarry Smith 13164b9ad928SBarry Smith Level: advanced 13174b9ad928SBarry Smith 131897177400SBarry Smith .seealso: PCMGSetLevels() 13194b9ad928SBarry Smith @*/ 13207087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 13214b9ad928SBarry Smith { 1322e952c529SMatthew G. Knepley PetscErrorCode ierr; 13234b9ad928SBarry Smith 13244b9ad928SBarry Smith PetscFunctionBegin; 13250700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13264482741eSBarry Smith PetscValidIntPointer(levels,2); 132797d33e41SMatthew G. Knepley *levels = 0; 132837b5128cSMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr); 13294b9ad928SBarry Smith PetscFunctionReturn(0); 13304b9ad928SBarry Smith } 13314b9ad928SBarry Smith 13324b9ad928SBarry Smith /*@ 133397177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 13344b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 13354b9ad928SBarry Smith 1336ad4df100SBarry Smith Logically Collective on PC 13374b9ad928SBarry Smith 13384b9ad928SBarry Smith Input Parameters: 13394b9ad928SBarry Smith + pc - the preconditioner context 13409dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 13419dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 13424b9ad928SBarry Smith 13434b9ad928SBarry Smith Options Database Key: 13444b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 13454b9ad928SBarry Smith additive, full, kaskade 13464b9ad928SBarry Smith 13474b9ad928SBarry Smith Level: advanced 13484b9ad928SBarry Smith 134997177400SBarry Smith .seealso: PCMGSetLevels() 13504b9ad928SBarry Smith @*/ 13517087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 13524b9ad928SBarry Smith { 1353f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 13544b9ad928SBarry Smith 13554b9ad928SBarry Smith PetscFunctionBegin; 13560700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1357c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 135831567311SBarry Smith mg->am = form; 13599dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 1360c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 1361c60c7ad4SBarry Smith PetscFunctionReturn(0); 1362c60c7ad4SBarry Smith } 1363c60c7ad4SBarry Smith 1364c60c7ad4SBarry Smith /*@ 1365c60c7ad4SBarry Smith PCMGGetType - Determines the form of multigrid to use: 1366c60c7ad4SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 1367c60c7ad4SBarry Smith 1368c60c7ad4SBarry Smith Logically Collective on PC 1369c60c7ad4SBarry Smith 1370c60c7ad4SBarry Smith Input Parameter: 1371c60c7ad4SBarry Smith . pc - the preconditioner context 1372c60c7ad4SBarry Smith 1373c60c7ad4SBarry Smith Output Parameter: 1374c60c7ad4SBarry Smith . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 1375c60c7ad4SBarry Smith 1376c60c7ad4SBarry Smith 1377c60c7ad4SBarry Smith Level: advanced 1378c60c7ad4SBarry Smith 1379c60c7ad4SBarry Smith .seealso: PCMGSetLevels() 1380c60c7ad4SBarry Smith @*/ 1381c60c7ad4SBarry Smith PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 1382c60c7ad4SBarry Smith { 1383c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1384c60c7ad4SBarry Smith 1385c60c7ad4SBarry Smith PetscFunctionBegin; 1386c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1387c60c7ad4SBarry Smith *type = mg->am; 13884b9ad928SBarry Smith PetscFunctionReturn(0); 13894b9ad928SBarry Smith } 13904b9ad928SBarry Smith 13914b9ad928SBarry Smith /*@ 13920d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 13934b9ad928SBarry Smith complicated cycling. 13944b9ad928SBarry Smith 1395ad4df100SBarry Smith Logically Collective on PC 13964b9ad928SBarry Smith 13974b9ad928SBarry Smith Input Parameters: 1398c2be2410SBarry Smith + pc - the multigrid context 1399c1cbb1deSBarry Smith - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 14004b9ad928SBarry Smith 14014b9ad928SBarry Smith Options Database Key: 1402c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired 14034b9ad928SBarry Smith 14044b9ad928SBarry Smith Level: advanced 14054b9ad928SBarry Smith 14060d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 14074b9ad928SBarry Smith @*/ 14087087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 14094b9ad928SBarry Smith { 1410f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1411f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 141279416396SBarry Smith PetscInt i,levels; 14134b9ad928SBarry Smith 14144b9ad928SBarry Smith PetscFunctionBegin; 14150700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 141669fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc,n,2); 14171a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1418f3fbd535SBarry Smith levels = mglevels[0]->levels; 14192fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 14204b9ad928SBarry Smith PetscFunctionReturn(0); 14214b9ad928SBarry Smith } 14224b9ad928SBarry Smith 14238cc2d5dfSBarry Smith /*@ 14248cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 14258cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 14268cc2d5dfSBarry Smith 1427ad4df100SBarry Smith Logically Collective on PC 14288cc2d5dfSBarry Smith 14298cc2d5dfSBarry Smith Input Parameters: 14308cc2d5dfSBarry Smith + pc - the multigrid context 14318cc2d5dfSBarry Smith - n - number of cycles (default is 1) 14328cc2d5dfSBarry Smith 14338cc2d5dfSBarry Smith Options Database Key: 1434e1bc860dSBarry Smith . -pc_mg_multiplicative_cycles n 14358cc2d5dfSBarry Smith 14368cc2d5dfSBarry Smith Level: advanced 14378cc2d5dfSBarry Smith 143895452b02SPatrick Sanan Notes: 143995452b02SPatrick Sanan This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 14408cc2d5dfSBarry Smith 14418cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 14428cc2d5dfSBarry Smith @*/ 14437087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 14448cc2d5dfSBarry Smith { 1445f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 14468cc2d5dfSBarry Smith 14478cc2d5dfSBarry Smith PetscFunctionBegin; 14480700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1449c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 14502a21e185SBarry Smith mg->cyclesperpcapply = n; 14518cc2d5dfSBarry Smith PetscFunctionReturn(0); 14528cc2d5dfSBarry Smith } 14538cc2d5dfSBarry Smith 14542134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1455fb15c04eSBarry Smith { 1456fb15c04eSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1457fb15c04eSBarry Smith 1458fb15c04eSBarry Smith PetscFunctionBegin; 14592134b1e4SBarry Smith mg->galerkin = use; 1460fb15c04eSBarry Smith PetscFunctionReturn(0); 1461fb15c04eSBarry Smith } 1462fb15c04eSBarry Smith 1463c2be2410SBarry Smith /*@ 146497177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 146582b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1466c2be2410SBarry Smith 1467ad4df100SBarry Smith Logically Collective on PC 1468c2be2410SBarry Smith 1469c2be2410SBarry Smith Input Parameters: 1470c91913d3SJed Brown + pc - the multigrid context 14712134b1e4SBarry Smith - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1472c2be2410SBarry Smith 1473c2be2410SBarry Smith Options Database Key: 14742134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> 1475c2be2410SBarry Smith 1476c2be2410SBarry Smith Level: intermediate 1477c2be2410SBarry Smith 147895452b02SPatrick Sanan Notes: 147995452b02SPatrick Sanan Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 14801c1aac46SBarry Smith use the PCMG construction of the coarser grids. 14811c1aac46SBarry Smith 14822134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType 14833fc8bf9cSBarry Smith 1484c2be2410SBarry Smith @*/ 14852134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1486c2be2410SBarry Smith { 1487fb15c04eSBarry Smith PetscErrorCode ierr; 1488c2be2410SBarry Smith 1489c2be2410SBarry Smith PetscFunctionBegin; 14900700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14912134b1e4SBarry Smith ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr); 1492c2be2410SBarry Smith PetscFunctionReturn(0); 1493c2be2410SBarry Smith } 1494c2be2410SBarry Smith 14953fc8bf9cSBarry Smith /*@ 14963fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 149782b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 14983fc8bf9cSBarry Smith 14993fc8bf9cSBarry Smith Not Collective 15003fc8bf9cSBarry Smith 15013fc8bf9cSBarry Smith Input Parameter: 15023fc8bf9cSBarry Smith . pc - the multigrid context 15033fc8bf9cSBarry Smith 15043fc8bf9cSBarry Smith Output Parameter: 15052134b1e4SBarry 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 15063fc8bf9cSBarry Smith 15073fc8bf9cSBarry Smith Level: intermediate 15083fc8bf9cSBarry Smith 15092134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType 15103fc8bf9cSBarry Smith 15113fc8bf9cSBarry Smith @*/ 15122134b1e4SBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 15133fc8bf9cSBarry Smith { 1514f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 15153fc8bf9cSBarry Smith 15163fc8bf9cSBarry Smith PetscFunctionBegin; 15170700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 15182134b1e4SBarry Smith *galerkin = mg->galerkin; 15193fc8bf9cSBarry Smith PetscFunctionReturn(0); 15203fc8bf9cSBarry Smith } 15213fc8bf9cSBarry Smith 1522f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt) 1523f3b08a26SMatthew G. Knepley { 1524f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1525f3b08a26SMatthew G. Knepley 1526f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1527f3b08a26SMatthew G. Knepley mg->adaptInterpolation = adapt; 1528f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1529f3b08a26SMatthew G. Knepley } 1530f3b08a26SMatthew G. Knepley 1531f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt) 1532f3b08a26SMatthew G. Knepley { 1533f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1534f3b08a26SMatthew G. Knepley 1535f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1536f3b08a26SMatthew G. Knepley *adapt = mg->adaptInterpolation; 1537f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1538f3b08a26SMatthew G. Knepley } 1539f3b08a26SMatthew G. Knepley 154041b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr) 154141b6fd38SMatthew G. Knepley { 154241b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 154341b6fd38SMatthew G. Knepley 154441b6fd38SMatthew G. Knepley PetscFunctionBegin; 154541b6fd38SMatthew G. Knepley mg->compatibleRelaxation = cr; 154641b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 154741b6fd38SMatthew G. Knepley } 154841b6fd38SMatthew G. Knepley 154941b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr) 155041b6fd38SMatthew G. Knepley { 155141b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 155241b6fd38SMatthew G. Knepley 155341b6fd38SMatthew G. Knepley PetscFunctionBegin; 155441b6fd38SMatthew G. Knepley *cr = mg->compatibleRelaxation; 155541b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 155641b6fd38SMatthew G. Knepley } 155741b6fd38SMatthew G. Knepley 1558f3b08a26SMatthew G. Knepley /*@ 1559f3b08a26SMatthew 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. 1560f3b08a26SMatthew G. Knepley 1561f3b08a26SMatthew G. Knepley Logically Collective on PC 1562f3b08a26SMatthew G. Knepley 1563f3b08a26SMatthew G. Knepley Input Parameters: 1564f3b08a26SMatthew G. Knepley + pc - the multigrid context 1565f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator 1566f3b08a26SMatthew G. Knepley 1567f3b08a26SMatthew G. Knepley Options Database Keys: 1568f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp - Turn on adaptation 1569f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension 1570f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector 1571f3b08a26SMatthew G. Knepley 1572f3b08a26SMatthew G. Knepley Level: intermediate 1573f3b08a26SMatthew G. Knepley 1574f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1575f3b08a26SMatthew G. Knepley .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin() 1576f3b08a26SMatthew G. Knepley @*/ 1577f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt) 1578f3b08a26SMatthew G. Knepley { 1579f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1580f3b08a26SMatthew G. Knepley 1581f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1582f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1583f3b08a26SMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));CHKERRQ(ierr); 1584f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1585f3b08a26SMatthew G. Knepley } 1586f3b08a26SMatthew G. Knepley 1587f3b08a26SMatthew G. Knepley /*@ 1588f3b08a26SMatthew 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. 1589f3b08a26SMatthew G. Knepley 1590f3b08a26SMatthew G. Knepley Not collective 1591f3b08a26SMatthew G. Knepley 1592f3b08a26SMatthew G. Knepley Input Parameter: 1593f3b08a26SMatthew G. Knepley . pc - the multigrid context 1594f3b08a26SMatthew G. Knepley 1595f3b08a26SMatthew G. Knepley Output Parameter: 1596f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator 1597f3b08a26SMatthew G. Knepley 1598f3b08a26SMatthew G. Knepley Level: intermediate 1599f3b08a26SMatthew G. Knepley 1600f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1601f3b08a26SMatthew G. Knepley .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin() 1602f3b08a26SMatthew G. Knepley @*/ 1603f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt) 1604f3b08a26SMatthew G. Knepley { 1605f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1606f3b08a26SMatthew G. Knepley 1607f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1608f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1609f3b08a26SMatthew G. Knepley PetscValidPointer(adapt, 2); 1610f3b08a26SMatthew G. Knepley ierr = PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));CHKERRQ(ierr); 1611f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1612f3b08a26SMatthew G. Knepley } 1613f3b08a26SMatthew G. Knepley 16144b9ad928SBarry Smith /*@ 161541b6fd38SMatthew G. Knepley PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation. 161641b6fd38SMatthew G. Knepley 161741b6fd38SMatthew G. Knepley Logically Collective on PC 161841b6fd38SMatthew G. Knepley 161941b6fd38SMatthew G. Knepley Input Parameters: 162041b6fd38SMatthew G. Knepley + pc - the multigrid context 162141b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation 162241b6fd38SMatthew G. Knepley 162341b6fd38SMatthew G. Knepley Options Database Keys: 162441b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation 162541b6fd38SMatthew G. Knepley 162641b6fd38SMatthew G. Knepley Level: intermediate 162741b6fd38SMatthew G. Knepley 162841b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin 162941b6fd38SMatthew G. Knepley .seealso: PCMGGetAdaptCR(), PCMGSetAdaptInterpolation(), PCMGSetGalerkin() 163041b6fd38SMatthew G. Knepley @*/ 163141b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr) 163241b6fd38SMatthew G. Knepley { 163341b6fd38SMatthew G. Knepley PetscErrorCode ierr; 163441b6fd38SMatthew G. Knepley 163541b6fd38SMatthew G. Knepley PetscFunctionBegin; 163641b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 163741b6fd38SMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr));CHKERRQ(ierr); 163841b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 163941b6fd38SMatthew G. Knepley } 164041b6fd38SMatthew G. Knepley 164141b6fd38SMatthew G. Knepley /*@ 164241b6fd38SMatthew G. Knepley PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation. 164341b6fd38SMatthew G. Knepley 164441b6fd38SMatthew G. Knepley Not collective 164541b6fd38SMatthew G. Knepley 164641b6fd38SMatthew G. Knepley Input Parameter: 164741b6fd38SMatthew G. Knepley . pc - the multigrid context 164841b6fd38SMatthew G. Knepley 164941b6fd38SMatthew G. Knepley Output Parameter: 165041b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion 165141b6fd38SMatthew G. Knepley 165241b6fd38SMatthew G. Knepley Level: intermediate 165341b6fd38SMatthew G. Knepley 165441b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin 165541b6fd38SMatthew G. Knepley .seealso: PCMGSetAdaptCR(), PCMGGetAdaptInterpolation(), PCMGSetGalerkin() 165641b6fd38SMatthew G. Knepley @*/ 165741b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr) 165841b6fd38SMatthew G. Knepley { 165941b6fd38SMatthew G. Knepley PetscErrorCode ierr; 166041b6fd38SMatthew G. Knepley 166141b6fd38SMatthew G. Knepley PetscFunctionBegin; 166241b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 166341b6fd38SMatthew G. Knepley PetscValidPointer(cr, 2); 166441b6fd38SMatthew G. Knepley ierr = PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr));CHKERRQ(ierr); 166541b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 166641b6fd38SMatthew G. Knepley } 166741b6fd38SMatthew G. Knepley 166841b6fd38SMatthew G. Knepley /*@ 166906792cafSBarry Smith PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1670710315b6SLawrence Mitchell on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1671710315b6SLawrence Mitchell pre- and post-smoothing steps. 167206792cafSBarry Smith 167306792cafSBarry Smith Logically Collective on PC 167406792cafSBarry Smith 167506792cafSBarry Smith Input Parameters: 167606792cafSBarry Smith + mg - the multigrid context 167706792cafSBarry Smith - n - the number of smoothing steps 167806792cafSBarry Smith 167906792cafSBarry Smith Options Database Key: 1680a2b725a8SWilliam Gropp . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 168106792cafSBarry Smith 168206792cafSBarry Smith Level: advanced 168306792cafSBarry Smith 168495452b02SPatrick Sanan Notes: 168595452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 168606792cafSBarry Smith there is no separate smooth up on the coarsest grid. 168706792cafSBarry Smith 1688710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp() 168906792cafSBarry Smith @*/ 169006792cafSBarry Smith PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 169106792cafSBarry Smith { 169206792cafSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 169306792cafSBarry Smith PC_MG_Levels **mglevels = mg->levels; 169406792cafSBarry Smith PetscErrorCode ierr; 169506792cafSBarry Smith PetscInt i,levels; 169606792cafSBarry Smith 169706792cafSBarry Smith PetscFunctionBegin; 169806792cafSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 169906792cafSBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 17001a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 170106792cafSBarry Smith levels = mglevels[0]->levels; 170206792cafSBarry Smith 170306792cafSBarry Smith for (i=1; i<levels; i++) { 170406792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 170506792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 170606792cafSBarry Smith mg->default_smoothu = n; 170706792cafSBarry Smith mg->default_smoothd = n; 170806792cafSBarry Smith } 170906792cafSBarry Smith PetscFunctionReturn(0); 171006792cafSBarry Smith } 171106792cafSBarry Smith 1712f442ab6aSBarry Smith /*@ 17138e5aa403SBarry Smith PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels 1714710315b6SLawrence Mitchell and adds the suffix _up to the options name 1715f442ab6aSBarry Smith 1716f442ab6aSBarry Smith Logically Collective on PC 1717f442ab6aSBarry Smith 1718f442ab6aSBarry Smith Input Parameters: 1719f442ab6aSBarry Smith . pc - the preconditioner context 1720f442ab6aSBarry Smith 1721f442ab6aSBarry Smith Options Database Key: 1722f442ab6aSBarry Smith . -pc_mg_distinct_smoothup 1723f442ab6aSBarry Smith 1724f442ab6aSBarry Smith Level: advanced 1725f442ab6aSBarry Smith 172695452b02SPatrick Sanan Notes: 172795452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 1728f442ab6aSBarry Smith there is no separate smooth up on the coarsest grid. 1729f442ab6aSBarry Smith 1730710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth() 1731f442ab6aSBarry Smith @*/ 1732f442ab6aSBarry Smith PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1733f442ab6aSBarry Smith { 1734f442ab6aSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1735f442ab6aSBarry Smith PC_MG_Levels **mglevels = mg->levels; 1736f442ab6aSBarry Smith PetscErrorCode ierr; 1737f442ab6aSBarry Smith PetscInt i,levels; 1738f442ab6aSBarry Smith KSP subksp; 1739f442ab6aSBarry Smith 1740f442ab6aSBarry Smith PetscFunctionBegin; 1741f442ab6aSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1742f442ab6aSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1743f442ab6aSBarry Smith levels = mglevels[0]->levels; 1744f442ab6aSBarry Smith 1745f442ab6aSBarry Smith for (i=1; i<levels; i++) { 1746710315b6SLawrence Mitchell const char *prefix = NULL; 1747f442ab6aSBarry Smith /* make sure smoother up and down are different */ 1748f442ab6aSBarry Smith ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr); 1749710315b6SLawrence Mitchell ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr); 1750710315b6SLawrence Mitchell ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr); 1751f442ab6aSBarry Smith ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr); 1752f442ab6aSBarry Smith } 1753f442ab6aSBarry Smith PetscFunctionReturn(0); 1754f442ab6aSBarry Smith } 1755f442ab6aSBarry Smith 175607a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 175707a4832bSFande Kong PetscErrorCode PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[]) 175807a4832bSFande Kong { 175907a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 176007a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 176107a4832bSFande Kong Mat *mat; 176207a4832bSFande Kong PetscInt l; 176307a4832bSFande Kong PetscErrorCode ierr; 176407a4832bSFande Kong 176507a4832bSFande Kong PetscFunctionBegin; 176607a4832bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 176707a4832bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 176807a4832bSFande Kong for (l=1; l< mg->nlevels; l++) { 176907a4832bSFande Kong mat[l-1] = mglevels[l]->interpolate; 177007a4832bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr); 177107a4832bSFande Kong } 177207a4832bSFande Kong *num_levels = mg->nlevels; 177307a4832bSFande Kong *interpolations = mat; 177407a4832bSFande Kong PetscFunctionReturn(0); 177507a4832bSFande Kong } 177607a4832bSFande Kong 177707a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 177807a4832bSFande Kong PetscErrorCode PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[]) 177907a4832bSFande Kong { 178007a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 178107a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 178207a4832bSFande Kong PetscInt l; 178307a4832bSFande Kong Mat *mat; 178407a4832bSFande Kong PetscErrorCode ierr; 178507a4832bSFande Kong 178607a4832bSFande Kong PetscFunctionBegin; 178707a4832bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 178807a4832bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 178907a4832bSFande Kong for (l=0; l<mg->nlevels-1; l++) { 179007a4832bSFande Kong ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr); 179107a4832bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr); 179207a4832bSFande Kong } 179307a4832bSFande Kong *num_levels = mg->nlevels; 179407a4832bSFande Kong *coarseOperators = mat; 179507a4832bSFande Kong PetscFunctionReturn(0); 179607a4832bSFande Kong } 179707a4832bSFande Kong 1798f3b08a26SMatthew G. Knepley /*@C 1799f3b08a26SMatthew G. Knepley PCMGRegisterCoarseSpaceConstructor - Adds a method to the PCMG package for coarse space construction. 1800f3b08a26SMatthew G. Knepley 1801f3b08a26SMatthew G. Knepley Not collective 1802f3b08a26SMatthew G. Knepley 1803f3b08a26SMatthew G. Knepley Input Parameters: 1804f3b08a26SMatthew G. Knepley + name - name of the constructor 1805f3b08a26SMatthew G. Knepley - function - constructor routine 1806f3b08a26SMatthew G. Knepley 1807f3b08a26SMatthew G. Knepley Notes: 1808f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1809f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp) 1810f3b08a26SMatthew G. Knepley $ pc - The PC object 1811f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1812f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1813f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1814f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1815f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1816f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1817f3b08a26SMatthew G. Knepley 1818f3b08a26SMatthew G. Knepley Level: advanced 1819f3b08a26SMatthew G. Knepley 1820f3b08a26SMatthew G. Knepley .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister() 1821f3b08a26SMatthew G. Knepley @*/ 1822f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)) 1823f3b08a26SMatthew G. Knepley { 1824f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1825f3b08a26SMatthew G. Knepley 1826f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1827f3b08a26SMatthew G. Knepley ierr = PCInitializePackage();CHKERRQ(ierr); 1828f3b08a26SMatthew G. Knepley ierr = PetscFunctionListAdd(&PCMGCoarseList,name,function);CHKERRQ(ierr); 1829f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1830f3b08a26SMatthew G. Knepley } 1831f3b08a26SMatthew G. Knepley 1832f3b08a26SMatthew G. Knepley /*@C 1833f3b08a26SMatthew G. Knepley PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method. 1834f3b08a26SMatthew G. Knepley 1835f3b08a26SMatthew G. Knepley Not collective 1836f3b08a26SMatthew G. Knepley 1837f3b08a26SMatthew G. Knepley Input Parameter: 1838f3b08a26SMatthew G. Knepley . name - name of the constructor 1839f3b08a26SMatthew G. Knepley 1840f3b08a26SMatthew G. Knepley Output Parameter: 1841f3b08a26SMatthew G. Knepley . function - constructor routine 1842f3b08a26SMatthew G. Knepley 1843f3b08a26SMatthew G. Knepley Notes: 1844f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1845f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp) 1846f3b08a26SMatthew G. Knepley $ pc - The PC object 1847f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1848f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1849f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1850f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1851f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1852f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1853f3b08a26SMatthew G. Knepley 1854f3b08a26SMatthew G. Knepley Level: advanced 1855f3b08a26SMatthew G. Knepley 1856f3b08a26SMatthew G. Knepley .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister() 1857f3b08a26SMatthew G. Knepley @*/ 1858f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)) 1859f3b08a26SMatthew G. Knepley { 1860f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1861f3b08a26SMatthew G. Knepley 1862f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1863f3b08a26SMatthew G. Knepley ierr = PetscFunctionListFind(PCMGCoarseList,name,function);CHKERRQ(ierr); 1864f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1865f3b08a26SMatthew G. Knepley } 1866f3b08a26SMatthew G. Knepley 18674b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 18684b9ad928SBarry Smith 18693b09bd56SBarry Smith /*MC 1870ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 18713b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 18723b09bd56SBarry Smith 18733b09bd56SBarry Smith Options Database Keys: 18743b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1875391689abSStefano Zampini . -pc_mg_cycle_type <v,w> - provide the cycle desired 18768c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 18773b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 1878710315b6SLawrence Mitchell . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 18792134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 18808cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 18818cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1882e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 18838cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 18848cf6037eSBarry Smith to the binary output file called binaryoutput 18853b09bd56SBarry Smith 188695452b02SPatrick Sanan Notes: 18870b3c753eSRichard 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 18883b09bd56SBarry Smith 18898cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 18908cf6037eSBarry Smith 189123067569SBarry Smith When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 189223067569SBarry Smith is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 189323067569SBarry Smith (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 189423067569SBarry Smith residual is computed at the end of each cycle. 189523067569SBarry Smith 18963b09bd56SBarry Smith Level: intermediate 18973b09bd56SBarry Smith 18988cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1899710315b6SLawrence Mitchell PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), 1900710315b6SLawrence Mitchell PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 190197177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 19020d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 19033b09bd56SBarry Smith M*/ 19043b09bd56SBarry Smith 19058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 19064b9ad928SBarry Smith { 1907f3fbd535SBarry Smith PC_MG *mg; 1908f3fbd535SBarry Smith PetscErrorCode ierr; 1909f3fbd535SBarry Smith 19104b9ad928SBarry Smith PetscFunctionBegin; 1911b00a9115SJed Brown ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1912f3fbd535SBarry Smith pc->data = (void*)mg; 1913f3fbd535SBarry Smith mg->nlevels = -1; 191410eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 19152134b1e4SBarry Smith mg->galerkin = PC_MG_GALERKIN_NONE; 1916f3b08a26SMatthew G. Knepley mg->adaptInterpolation = PETSC_FALSE; 1917f3b08a26SMatthew G. Knepley mg->Nc = -1; 1918f3b08a26SMatthew G. Knepley mg->eigenvalue = -1; 1919f3fbd535SBarry Smith 192037a44384SMark Adams pc->useAmat = PETSC_TRUE; 192137a44384SMark Adams 19224b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 1923fcb023d4SJed Brown pc->ops->applytranspose = PCApplyTranspose_MG; 192430b0564aSStefano Zampini pc->ops->matapply = PCMatApply_MG; 19254b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1926a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 19274b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 19284b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 19294b9ad928SBarry Smith pc->ops->view = PCView_MG; 1930fb15c04eSBarry Smith 1931f3b08a26SMatthew G. Knepley ierr = PetscObjectComposedDataRegister(&mg->eigenvalue);CHKERRQ(ierr); 1932fb15c04eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 193397d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr); 193497d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr); 1935fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);CHKERRQ(ierr); 1936fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);CHKERRQ(ierr); 1937f3b08a26SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);CHKERRQ(ierr); 1938f3b08a26SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);CHKERRQ(ierr); 193941b6fd38SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG);CHKERRQ(ierr); 194041b6fd38SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG);CHKERRQ(ierr); 19414b9ad928SBarry Smith PetscFunctionReturn(0); 19424b9ad928SBarry Smith } 1943