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 15*30b0564aSStefano 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) { 25*30b0564aSStefano Zampini if (matapp) { 26*30b0564aSStefano Zampini ierr = KSPMatSolve(mglevels->smoothd,mglevels->B,mglevels->X);CHKERRQ(ierr); /* pre-smooth */ 27*30b0564aSStefano Zampini ierr = KSPCheckSolve(mglevels->smoothd,pc,NULL);CHKERRQ(ierr); 28*30b0564aSStefano 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); 31*30b0564aSStefano Zampini } 32fcb023d4SJed Brown } else { 33*30b0564aSStefano 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);} 40*30b0564aSStefano Zampini if (matapp && !mglevels->R) { 41*30b0564aSStefano Zampini ierr = MatDuplicate(mglevels->B,MAT_DO_NOT_COPY_VALUES,&mglevels->R);CHKERRQ(ierr); 42*30b0564aSStefano Zampini } 43fcb023d4SJed Brown if (!transpose) { 44*30b0564aSStefano Zampini if (matapp) { ierr = (*mglevels->matresidual)(mglevels->A,mglevels->B,mglevels->X,mglevels->R);CHKERRQ(ierr); } 45*30b0564aSStefano Zampini else { ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); } 46fcb023d4SJed Brown } else { 47*30b0564aSStefano Zampini if (matapp) { ierr = (*mglevels->matresidualtranspose)(mglevels->A,mglevels->B,mglevels->X,mglevels->R);CHKERRQ(ierr); } 48*30b0564aSStefano 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) { 71*30b0564aSStefano Zampini if (matapp) { ierr = MatMatRestrict(mglevels->restrct,mglevels->R,&mgc->B);CHKERRQ(ierr); } 72*30b0564aSStefano Zampini else { ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); } 73fcb023d4SJed Brown } else { 74*30b0564aSStefano Zampini if (matapp) { ierr = MatMatRestrict(mglevels->interpolate,mglevels->R,&mgc->B);CHKERRQ(ierr); } 75*30b0564aSStefano 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);} 78*30b0564aSStefano Zampini if (matapp) { 79*30b0564aSStefano Zampini if (!mgc->X) { 80*30b0564aSStefano Zampini ierr = MatDuplicate(mgc->B,MAT_DO_NOT_COPY_VALUES,&mgc->X);CHKERRQ(ierr); 81*30b0564aSStefano Zampini } else { 82*30b0564aSStefano Zampini ierr = MatZeroEntries(mgc->X);CHKERRQ(ierr); 83*30b0564aSStefano Zampini } 84*30b0564aSStefano Zampini } else { 85*30b0564aSStefano Zampini ierr = VecZeroEntries(mgc->x);CHKERRQ(ierr); 86*30b0564aSStefano Zampini } 874b9ad928SBarry Smith while (cycles--) { 88*30b0564aSStefano 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) { 92*30b0564aSStefano Zampini if (matapp) { ierr = MatMatInterpolateAdd(mglevels->interpolate,mgc->X,mglevels->X,&mglevels->X);CHKERRQ(ierr); } 93*30b0564aSStefano 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) { 100*30b0564aSStefano Zampini if (matapp) { 101*30b0564aSStefano Zampini ierr = KSPMatSolve(mglevels->smoothu,mglevels->B,mglevels->X);CHKERRQ(ierr); /* post smooth */ 102*30b0564aSStefano Zampini ierr = KSPCheckSolve(mglevels->smoothu,pc,NULL);CHKERRQ(ierr); 103*30b0564aSStefano 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); 106*30b0564aSStefano Zampini } 107fcb023d4SJed Brown } else { 108*30b0564aSStefano 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) { 113*30b0564aSStefano 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++) { 192*30b0564aSStefano 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); 215*30b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i+1]->R);CHKERRQ(ierr); 216*30b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i]->B);CHKERRQ(ierr); 217*30b0564aSStefano 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); 230*30b0564aSStefano 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 462d5a8f7e6SBarry Smith you must pass MPI_COMM_NULL. 46397d33e41SMatthew G. Knepley 46497d33e41SMatthew G. Knepley Level: intermediate 46597d33e41SMatthew G. Knepley 46697d33e41SMatthew G. Knepley Notes: 46797d33e41SMatthew G. Knepley If the number of levels is one then the multigrid uses the -mg_levels prefix 46897d33e41SMatthew G. Knepley for setting the level options rather than the -mg_coarse prefix. 46997d33e41SMatthew G. Knepley 470d5a8f7e6SBarry Smith You can free the information in comms after this routine is called. 471d5a8f7e6SBarry Smith 472d5a8f7e6SBarry Smith The array of MPI communicators must contain MPI_COMM_NULL for those ranks that at each level 473d5a8f7e6SBarry Smith are not participating in the coarser solve. For example, with 2 levels and 1 and 2 ranks on 474d5a8f7e6SBarry Smith the two levels, rank 0 in the original communicator will pass in an array of 2 communicators 475d5a8f7e6SBarry Smith of size 2 and 1, while rank 1 in the original communicator will pass in array of 2 communicators 476d5a8f7e6SBarry Smith the first of size 2 and the second of value MPI_COMM_NULL since the rank 1 does not participate 477d5a8f7e6SBarry Smith in the coarse grid solve. 478d5a8f7e6SBarry Smith 479d5a8f7e6SBarry Smith Since each coarser level may have a new MPI_Comm with fewer ranks than the previous, one 480d5a8f7e6SBarry Smith must take special care in providing the restriction and interpolation operation. We recommend 481d5a8f7e6SBarry Smith providing these as two step operations; first perform a standard restriction or interpolation on 482d5a8f7e6SBarry Smith the full number of ranks for that level and then use an MPI call to copy the resulting vector 483d5a8f7e6SBarry Smith array entries (after calls to VecGetArray()) to the smaller or larger number of ranks, not in both 484d5a8f7e6SBarry Smith cases the MPI calls must be made on the larger of the two communicators. Traditional MPI send and 485d5a8f7e6SBarry Smith recieves or MPI_AlltoAllv() could be used to do the reshuffling of the vector entries. 486d5a8f7e6SBarry Smith 487d5a8f7e6SBarry Smith 48897d33e41SMatthew G. Knepley .seealso: PCMGSetType(), PCMGGetLevels() 48997d33e41SMatthew G. Knepley @*/ 49097d33e41SMatthew G. Knepley PetscErrorCode PCMGSetLevels(PC pc,PetscInt levels,MPI_Comm *comms) 49197d33e41SMatthew G. Knepley { 49297d33e41SMatthew G. Knepley PetscErrorCode ierr; 49397d33e41SMatthew G. Knepley 49497d33e41SMatthew G. Knepley PetscFunctionBegin; 49597d33e41SMatthew G. Knepley PetscValidHeaderSpecific(pc,PC_CLASSID,1); 49697d33e41SMatthew G. Knepley if (comms) PetscValidPointer(comms,3); 49737b5128cSMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetLevels_C",(PC,PetscInt,MPI_Comm*),(pc,levels,comms));CHKERRQ(ierr); 49897d33e41SMatthew G. Knepley PetscFunctionReturn(0); 49997d33e41SMatthew G. Knepley } 50097d33e41SMatthew G. Knepley 501c07bf074SBarry Smith 502c07bf074SBarry Smith PetscErrorCode PCDestroy_MG(PC pc) 503c07bf074SBarry Smith { 504c07bf074SBarry Smith PetscErrorCode ierr; 505a06653b4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 506a06653b4SBarry Smith PC_MG_Levels **mglevels = mg->levels; 507a06653b4SBarry Smith PetscInt i,n; 508c07bf074SBarry Smith 509c07bf074SBarry Smith PetscFunctionBegin; 510a06653b4SBarry Smith ierr = PCReset_MG(pc);CHKERRQ(ierr); 511a06653b4SBarry Smith if (mglevels) { 512a06653b4SBarry Smith n = mglevels[0]->levels; 513a06653b4SBarry Smith for (i=0; i<n; i++) { 514a06653b4SBarry Smith if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 5156bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothd);CHKERRQ(ierr); 516a06653b4SBarry Smith } 5176bf464f9SBarry Smith ierr = KSPDestroy(&mglevels[i]->smoothu);CHKERRQ(ierr); 51841b6fd38SMatthew G. Knepley ierr = KSPDestroy(&mglevels[i]->cr);CHKERRQ(ierr); 519a06653b4SBarry Smith ierr = PetscFree(mglevels[i]);CHKERRQ(ierr); 520a06653b4SBarry Smith } 521a06653b4SBarry Smith ierr = PetscFree(mg->levels);CHKERRQ(ierr); 522a06653b4SBarry Smith } 523c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 524fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",NULL);CHKERRQ(ierr); 525fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",NULL);CHKERRQ(ierr); 526f3fbd535SBarry Smith PetscFunctionReturn(0); 527f3fbd535SBarry Smith } 528f3fbd535SBarry Smith 529f3fbd535SBarry Smith 530f3fbd535SBarry Smith /* 531f3fbd535SBarry Smith PCApply_MG - Runs either an additive, multiplicative, Kaskadic 532f3fbd535SBarry Smith or full cycle of multigrid. 533f3fbd535SBarry Smith 534f3fbd535SBarry Smith Note: 535f3fbd535SBarry Smith A simple wrapper which calls PCMGMCycle(),PCMGACycle(), or PCMGFCycle(). 536f3fbd535SBarry Smith */ 537*30b0564aSStefano Zampini static PetscErrorCode PCApply_MG_Internal(PC pc,Vec b,Vec x,Mat B,Mat X,PetscBool transpose) 538f3fbd535SBarry Smith { 539f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 540f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 541f3fbd535SBarry Smith PetscErrorCode ierr; 542391689abSStefano Zampini PC tpc; 543f3fbd535SBarry Smith PetscInt levels = mglevels[0]->levels,i; 544*30b0564aSStefano Zampini PetscBool changeu,changed,matapp; 545f3fbd535SBarry Smith 546f3fbd535SBarry Smith PetscFunctionBegin; 547*30b0564aSStefano Zampini matapp = (PetscBool)(B && X); 548a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} 549e1d8e5deSBarry Smith /* When the DM is supplying the matrix then it will not exist until here */ 550298cc208SBarry Smith for (i=0; i<levels; i++) { 551e1d8e5deSBarry Smith if (!mglevels[i]->A) { 55223ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); 553298cc208SBarry Smith ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); 554e1d8e5deSBarry Smith } 555e1d8e5deSBarry Smith } 556e1d8e5deSBarry Smith 557391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothd,&tpc);CHKERRQ(ierr); 558391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changed);CHKERRQ(ierr); 559391689abSStefano Zampini ierr = KSPGetPC(mglevels[levels-1]->smoothu,&tpc);CHKERRQ(ierr); 560391689abSStefano Zampini ierr = PCPreSolveChangeRHS(tpc,&changeu);CHKERRQ(ierr); 561391689abSStefano Zampini if (!changeu && !changed) { 562*30b0564aSStefano Zampini if (matapp) { 563*30b0564aSStefano Zampini ierr = MatDestroy(&mglevels[levels-1]->B);CHKERRQ(ierr); 564*30b0564aSStefano Zampini mglevels[levels-1]->B = B; 565*30b0564aSStefano Zampini } else { 566391689abSStefano Zampini ierr = VecDestroy(&mglevels[levels-1]->b);CHKERRQ(ierr); 567f3fbd535SBarry Smith mglevels[levels-1]->b = b; 568*30b0564aSStefano Zampini } 569391689abSStefano Zampini } else { /* if the smoother changes the rhs during PreSolve, we cannot use the input vector */ 570*30b0564aSStefano Zampini if (matapp) { 571*30b0564aSStefano Zampini if (mglevels[levels-1]->B) { 572*30b0564aSStefano Zampini PetscInt N1,N2; 573*30b0564aSStefano Zampini PetscBool flg; 574*30b0564aSStefano Zampini 575*30b0564aSStefano Zampini ierr = MatGetSize(mglevels[levels-1]->B,NULL,&N1);CHKERRQ(ierr); 576*30b0564aSStefano Zampini ierr = MatGetSize(B,NULL,&N2);CHKERRQ(ierr); 577*30b0564aSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)mglevels[levels-1]->B,((PetscObject)B)->type_name,&flg);CHKERRQ(ierr); 578*30b0564aSStefano Zampini if (N1 != N2 || !flg) { 579*30b0564aSStefano Zampini ierr = MatDestroy(&mglevels[levels-1]->B);CHKERRQ(ierr); 580*30b0564aSStefano Zampini } 581*30b0564aSStefano Zampini } 582*30b0564aSStefano Zampini if (!mglevels[levels-1]->B) { 583*30b0564aSStefano Zampini ierr = MatDuplicate(B,MAT_COPY_VALUES,&mglevels[levels-1]->B);CHKERRQ(ierr); 584*30b0564aSStefano Zampini } else { 585*30b0564aSStefano Zampini ierr = MatCopy(B,mglevels[levels-1]->B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 586*30b0564aSStefano Zampini } 587*30b0564aSStefano Zampini } else { 588391689abSStefano Zampini if (!mglevels[levels-1]->b) { 589391689abSStefano Zampini Vec *vec; 590391689abSStefano Zampini 591391689abSStefano Zampini ierr = KSPCreateVecs(mglevels[levels-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 592391689abSStefano Zampini mglevels[levels-1]->b = *vec; 5937ae21d43SStefano Zampini ierr = PetscFree(vec);CHKERRQ(ierr); 594391689abSStefano Zampini } 595391689abSStefano Zampini ierr = VecCopy(b,mglevels[levels-1]->b);CHKERRQ(ierr); 596391689abSStefano Zampini } 597*30b0564aSStefano Zampini } 598*30b0564aSStefano Zampini if (matapp) { mglevels[levels-1]->X = X; } 599*30b0564aSStefano Zampini else { mglevels[levels-1]->x = x; } 600*30b0564aSStefano Zampini 601*30b0564aSStefano Zampini /* If coarser Xs are present, it means we have already block applied the PC at least once 602*30b0564aSStefano Zampini Reset operators if sizes/type do no match */ 603*30b0564aSStefano Zampini if (matapp && levels > 1 && mglevels[levels-2]->X) { 604*30b0564aSStefano Zampini PetscInt Xc,Bc; 605*30b0564aSStefano Zampini PetscBool flg; 606*30b0564aSStefano Zampini 607*30b0564aSStefano Zampini ierr = MatGetSize(mglevels[levels-2]->X,NULL,&Xc);CHKERRQ(ierr); 608*30b0564aSStefano Zampini ierr = MatGetSize(mglevels[levels-1]->B,NULL,&Bc);CHKERRQ(ierr); 609*30b0564aSStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)mglevels[levels-2]->X,((PetscObject)mglevels[levels-1]->X)->type_name,&flg);CHKERRQ(ierr); 610*30b0564aSStefano Zampini if (Xc != Bc || !flg) { 611*30b0564aSStefano Zampini ierr = MatDestroy(&mglevels[levels-1]->R);CHKERRQ(ierr); 612*30b0564aSStefano Zampini for (i=0;i<levels-1;i++) { 613*30b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i]->R);CHKERRQ(ierr); 614*30b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i]->B);CHKERRQ(ierr); 615*30b0564aSStefano Zampini ierr = MatDestroy(&mglevels[i]->X);CHKERRQ(ierr); 616*30b0564aSStefano Zampini } 617*30b0564aSStefano Zampini } 618*30b0564aSStefano Zampini } 619391689abSStefano Zampini 62031567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 621*30b0564aSStefano Zampini if (matapp) { ierr = MatZeroEntries(X);CHKERRQ(ierr); } 622*30b0564aSStefano Zampini else { ierr = VecZeroEntries(x);CHKERRQ(ierr); } 62331567311SBarry Smith for (i=0; i<mg->cyclesperpcapply; i++) { 624*30b0564aSStefano Zampini ierr = PCMGMCycle_Private(pc,mglevels+levels-1,transpose,matapp,NULL);CHKERRQ(ierr); 625f3fbd535SBarry Smith } 6262fa5cd67SKarl Rupp } else if (mg->am == PC_MG_ADDITIVE) { 627*30b0564aSStefano Zampini ierr = PCMGACycle_Private(pc,mglevels,transpose,matapp);CHKERRQ(ierr); 6282fa5cd67SKarl Rupp } else if (mg->am == PC_MG_KASKADE) { 629*30b0564aSStefano Zampini ierr = PCMGKCycle_Private(pc,mglevels,transpose,matapp);CHKERRQ(ierr); 6302fa5cd67SKarl Rupp } else { 631*30b0564aSStefano Zampini ierr = PCMGFCycle_Private(pc,mglevels,transpose,matapp);CHKERRQ(ierr); 632f3fbd535SBarry Smith } 633a9db87a2SMatthew G Knepley if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} 634*30b0564aSStefano Zampini if (!changeu && !changed) { 635*30b0564aSStefano Zampini if (matapp) { mglevels[levels-1]->B = NULL; } 636*30b0564aSStefano Zampini else { mglevels[levels-1]->b = NULL; } 637*30b0564aSStefano Zampini } 638f3fbd535SBarry Smith PetscFunctionReturn(0); 639f3fbd535SBarry Smith } 640f3fbd535SBarry Smith 641fcb023d4SJed Brown static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) 642fcb023d4SJed Brown { 643fcb023d4SJed Brown PetscErrorCode ierr; 644fcb023d4SJed Brown 645fcb023d4SJed Brown PetscFunctionBegin; 646*30b0564aSStefano Zampini ierr = PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_FALSE);CHKERRQ(ierr); 647fcb023d4SJed Brown PetscFunctionReturn(0); 648fcb023d4SJed Brown } 649fcb023d4SJed Brown 650fcb023d4SJed Brown static PetscErrorCode PCApplyTranspose_MG(PC pc,Vec b,Vec x) 651fcb023d4SJed Brown { 652fcb023d4SJed Brown PetscErrorCode ierr; 653fcb023d4SJed Brown 654fcb023d4SJed Brown PetscFunctionBegin; 655*30b0564aSStefano Zampini ierr = PCApply_MG_Internal(pc,b,x,NULL,NULL,PETSC_TRUE);CHKERRQ(ierr); 656*30b0564aSStefano Zampini PetscFunctionReturn(0); 657*30b0564aSStefano Zampini } 658*30b0564aSStefano Zampini 659*30b0564aSStefano Zampini static PetscErrorCode PCMatApply_MG(PC pc,Mat b,Mat x) 660*30b0564aSStefano Zampini { 661*30b0564aSStefano Zampini PetscErrorCode ierr; 662*30b0564aSStefano Zampini 663*30b0564aSStefano Zampini PetscFunctionBegin; 664*30b0564aSStefano Zampini ierr = PCApply_MG_Internal(pc,NULL,NULL,b,x,PETSC_FALSE);CHKERRQ(ierr); 665fcb023d4SJed Brown PetscFunctionReturn(0); 666fcb023d4SJed Brown } 667f3fbd535SBarry Smith 6684416b707SBarry Smith PetscErrorCode PCSetFromOptions_MG(PetscOptionItems *PetscOptionsObject,PC pc) 669f3fbd535SBarry Smith { 670f3fbd535SBarry Smith PetscErrorCode ierr; 671710315b6SLawrence Mitchell PetscInt levels,cycles; 672f3b08a26SMatthew G. Knepley PetscBool flg, flg2; 673f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 6743d3eaba7SBarry Smith PC_MG_Levels **mglevels; 675f3fbd535SBarry Smith PCMGType mgtype; 676f3fbd535SBarry Smith PCMGCycleType mgctype; 6772134b1e4SBarry Smith PCMGGalerkinType gtype; 678f3fbd535SBarry Smith 679f3fbd535SBarry Smith PetscFunctionBegin; 6801d743356SStefano Zampini levels = PetscMax(mg->nlevels,1); 681e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Multigrid options");CHKERRQ(ierr); 682f3fbd535SBarry Smith ierr = PetscOptionsInt("-pc_mg_levels","Number of Levels","PCMGSetLevels",levels,&levels,&flg);CHKERRQ(ierr); 6831a97d7b6SStefano Zampini if (!flg && !mg->levels && pc->dm) { 684eb3f98d2SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 685eb3f98d2SBarry Smith levels++; 6863aeaf226SBarry Smith mg->usedmfornumberoflevels = PETSC_TRUE; 687eb3f98d2SBarry Smith } 6880298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 6893aeaf226SBarry Smith mglevels = mg->levels; 6903aeaf226SBarry Smith 691f3fbd535SBarry Smith mgctype = (PCMGCycleType) mglevels[0]->cycles; 692f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_cycle_type","V cycle or for W-cycle","PCMGSetCycleType",PCMGCycleTypes,(PetscEnum)mgctype,(PetscEnum*)&mgctype,&flg);CHKERRQ(ierr); 693f3fbd535SBarry Smith if (flg) { 694f3fbd535SBarry Smith ierr = PCMGSetCycleType(pc,mgctype);CHKERRQ(ierr); 6952fa5cd67SKarl Rupp } 6962134b1e4SBarry Smith gtype = mg->galerkin; 6972134b1e4SBarry Smith ierr = PetscOptionsEnum("-pc_mg_galerkin","Use Galerkin process to compute coarser operators","PCMGSetGalerkin",PCMGGalerkinTypes,(PetscEnum)gtype,(PetscEnum*)>ype,&flg);CHKERRQ(ierr); 6982134b1e4SBarry Smith if (flg) { 6992134b1e4SBarry Smith ierr = PCMGSetGalerkin(pc,gtype);CHKERRQ(ierr); 700f3fbd535SBarry Smith } 701f3b08a26SMatthew G. Knepley flg2 = PETSC_FALSE; 702f3b08a26SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_adapt_interp","Adapt interpolation using some coarse space","PCMGSetAdaptInterpolation",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr); 703f3b08a26SMatthew G. Knepley if (flg) {ierr = PCMGSetAdaptInterpolation(pc, flg2);CHKERRQ(ierr);} 704f3b08a26SMatthew 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); 705f3b08a26SMatthew 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); 706f3b08a26SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_mesp_monitor","Monitor the multilevel eigensolver","PCMGSetAdaptInterpolation",PETSC_FALSE,&mg->mespMonitor,&flg);CHKERRQ(ierr); 70741b6fd38SMatthew G. Knepley flg2 = PETSC_FALSE; 70841b6fd38SMatthew G. Knepley ierr = PetscOptionsBool("-pc_mg_adapt_cr","Monitor coarse space quality using Compatible Relaxation (CR)","PCMGSetAdaptCR",PETSC_FALSE,&flg2,&flg);CHKERRQ(ierr); 70941b6fd38SMatthew G. Knepley if (flg) {ierr = PCMGSetAdaptCR(pc, flg2);CHKERRQ(ierr);} 710f56b1265SBarry Smith flg = PETSC_FALSE; 7118e5aa403SBarry Smith ierr = PetscOptionsBool("-pc_mg_distinct_smoothup","Create separate smoothup KSP and append the prefix _up","PCMGSetDistinctSmoothUp",PETSC_FALSE,&flg,NULL);CHKERRQ(ierr); 712f442ab6aSBarry Smith if (flg) { 713f442ab6aSBarry Smith ierr = PCMGSetDistinctSmoothUp(pc);CHKERRQ(ierr); 714f442ab6aSBarry Smith } 71531567311SBarry Smith mgtype = mg->am; 716f3fbd535SBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)mgtype,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 717f3fbd535SBarry Smith if (flg) { 718f3fbd535SBarry Smith ierr = PCMGSetType(pc,mgtype);CHKERRQ(ierr); 719f3fbd535SBarry Smith } 72031567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 7215363c1e0SLisandro Dalcin ierr = PetscOptionsInt("-pc_mg_multiplicative_cycles","Number of cycles for each preconditioner step","PCMGMultiplicativeSetCycles",mg->cyclesperpcapply,&cycles,&flg);CHKERRQ(ierr); 722f3fbd535SBarry Smith if (flg) { 723f3fbd535SBarry Smith ierr = PCMGMultiplicativeSetCycles(pc,cycles);CHKERRQ(ierr); 724f3fbd535SBarry Smith } 725f3fbd535SBarry Smith } 726f3fbd535SBarry Smith flg = PETSC_FALSE; 7270298fd71SBarry Smith ierr = PetscOptionsBool("-pc_mg_log","Log times for each multigrid level","None",flg,&flg,NULL);CHKERRQ(ierr); 728f3fbd535SBarry Smith if (flg) { 729f3fbd535SBarry Smith PetscInt i; 730f3fbd535SBarry Smith char eventname[128]; 7311a97d7b6SStefano Zampini 732f3fbd535SBarry Smith levels = mglevels[0]->levels; 733f3fbd535SBarry Smith for (i=0; i<levels; i++) { 734f3fbd535SBarry Smith sprintf(eventname,"MGSetup Level %d",(int)i); 73563e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsetup);CHKERRQ(ierr); 736f3fbd535SBarry Smith sprintf(eventname,"MGSmooth Level %d",(int)i); 73763e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventsmoothsolve);CHKERRQ(ierr); 738f3fbd535SBarry Smith if (i) { 739f3fbd535SBarry Smith sprintf(eventname,"MGResid Level %d",(int)i); 74063e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventresidual);CHKERRQ(ierr); 741f3fbd535SBarry Smith sprintf(eventname,"MGInterp Level %d",(int)i); 74263e6d426SJed Brown ierr = PetscLogEventRegister(eventname,((PetscObject)pc)->classid,&mglevels[i]->eventinterprestrict);CHKERRQ(ierr); 743f3fbd535SBarry Smith } 744f3fbd535SBarry Smith } 745eec35531SMatthew G Knepley 746ec60ca73SSatish Balay #if defined(PETSC_USE_LOG) 747eec35531SMatthew G Knepley { 748eec35531SMatthew G Knepley const char *sname = "MG Apply"; 749eec35531SMatthew G Knepley PetscStageLog stageLog; 750eec35531SMatthew G Knepley PetscInt st; 751eec35531SMatthew G Knepley 752eec35531SMatthew G Knepley ierr = PetscLogGetStageLog(&stageLog);CHKERRQ(ierr); 753eec35531SMatthew G Knepley for (st = 0; st < stageLog->numStages; ++st) { 754eec35531SMatthew G Knepley PetscBool same; 755eec35531SMatthew G Knepley 756eec35531SMatthew G Knepley ierr = PetscStrcmp(stageLog->stageInfo[st].name, sname, &same);CHKERRQ(ierr); 7572fa5cd67SKarl Rupp if (same) mg->stageApply = st; 758eec35531SMatthew G Knepley } 759eec35531SMatthew G Knepley if (!mg->stageApply) { 760eec35531SMatthew G Knepley ierr = PetscLogStageRegister(sname, &mg->stageApply);CHKERRQ(ierr); 761eec35531SMatthew G Knepley } 762eec35531SMatthew G Knepley } 763ec60ca73SSatish Balay #endif 764f3fbd535SBarry Smith } 765f3fbd535SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 766f3b08a26SMatthew G. Knepley /* Check option consistency */ 767f3b08a26SMatthew G. Knepley ierr = PCMGGetGalerkin(pc, >ype);CHKERRQ(ierr); 768f3b08a26SMatthew G. Knepley ierr = PCMGGetAdaptInterpolation(pc, &flg);CHKERRQ(ierr); 769f3b08a26SMatthew 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"); 770f3fbd535SBarry Smith PetscFunctionReturn(0); 771f3fbd535SBarry Smith } 772f3fbd535SBarry Smith 7730a545947SLisandro Dalcin const char *const PCMGTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","PCMGType","PC_MG",NULL}; 7740a545947SLisandro Dalcin const char *const PCMGCycleTypes[] = {"invalid","v","w","PCMGCycleType","PC_MG_CYCLE",NULL}; 7750a545947SLisandro Dalcin const char *const PCMGGalerkinTypes[] = {"both","pmat","mat","none","external","PCMGGalerkinType","PC_MG_GALERKIN",NULL}; 776f3b08a26SMatthew G. Knepley const char *const PCMGCoarseSpaceTypes[] = {"polynomial","harmonic","eigenvector","generalized_eigenvector","PCMGCoarseSpaceType","PCMG_POLYNOMIAL",NULL}; 777f3fbd535SBarry Smith 7789804daf3SBarry Smith #include <petscdraw.h> 779f3fbd535SBarry Smith PetscErrorCode PCView_MG(PC pc,PetscViewer viewer) 780f3fbd535SBarry Smith { 781f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 782f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 783f3fbd535SBarry Smith PetscErrorCode ierr; 784e3deeeafSJed Brown PetscInt levels = mglevels ? mglevels[0]->levels : 0,i; 785219b1045SBarry Smith PetscBool iascii,isbinary,isdraw; 786f3fbd535SBarry Smith 787f3fbd535SBarry Smith PetscFunctionBegin; 788251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 7895b0b0462SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 790219b1045SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 791f3fbd535SBarry Smith if (iascii) { 792e3deeeafSJed Brown const char *cyclename = levels ? (mglevels[0]->cycles == PC_MG_CYCLE_V ? "v" : "w") : "unknown"; 793efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type is %s, levels=%D cycles=%s\n", PCMGTypes[mg->am],levels,cyclename);CHKERRQ(ierr); 79431567311SBarry Smith if (mg->am == PC_MG_MULTIPLICATIVE) { 79531567311SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Cycles per PCApply=%d\n",mg->cyclesperpcapply);CHKERRQ(ierr); 796f3fbd535SBarry Smith } 7972134b1e4SBarry Smith if (mg->galerkin == PC_MG_GALERKIN_BOTH) { 798f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 7992134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_PMAT) { 8002134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for pmat\n");CHKERRQ(ierr); 8012134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_MAT) { 8022134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid matrices for mat\n");CHKERRQ(ierr); 8032134b1e4SBarry Smith } else if (mg->galerkin == PC_MG_GALERKIN_EXTERNAL) { 8042134b1e4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Using externally compute Galerkin coarse grid matrices\n");CHKERRQ(ierr); 8054f66f45eSBarry Smith } else { 8064f66f45eSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid matrices\n");CHKERRQ(ierr); 807f3fbd535SBarry Smith } 8085adeb434SBarry Smith if (mg->view){ 8095adeb434SBarry Smith ierr = (*mg->view)(pc,viewer);CHKERRQ(ierr); 8105adeb434SBarry Smith } 811f3fbd535SBarry Smith for (i=0; i<levels; i++) { 812f3fbd535SBarry Smith if (!i) { 81389cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level -------------------------------\n",i);CHKERRQ(ierr); 814f3fbd535SBarry Smith } else { 81589cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 816f3fbd535SBarry Smith } 817f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 818f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 819f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 820f3fbd535SBarry Smith if (i && mglevels[i]->smoothd == mglevels[i]->smoothu) { 821f3fbd535SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 822f3fbd535SBarry Smith } else if (i) { 82389cce641SBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 824f3fbd535SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 825f3fbd535SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 826f3fbd535SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 827f3fbd535SBarry Smith } 82841b6fd38SMatthew G. Knepley if (i && mglevels[i]->cr) { 82941b6fd38SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer,"CR solver on level %D -------------------------------\n",i);CHKERRQ(ierr); 83041b6fd38SMatthew G. Knepley ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 83141b6fd38SMatthew G. Knepley ierr = KSPView(mglevels[i]->cr,viewer);CHKERRQ(ierr); 83241b6fd38SMatthew G. Knepley ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 83341b6fd38SMatthew G. Knepley } 834f3fbd535SBarry Smith } 8355b0b0462SBarry Smith } else if (isbinary) { 8365b0b0462SBarry Smith for (i=levels-1; i>=0; i--) { 8375b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 8385b0b0462SBarry Smith if (i && mglevels[i]->smoothd != mglevels[i]->smoothu) { 8395b0b0462SBarry Smith ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 8405b0b0462SBarry Smith } 8415b0b0462SBarry Smith } 842219b1045SBarry Smith } else if (isdraw) { 843219b1045SBarry Smith PetscDraw draw; 844b4375e8dSPeter Brune PetscReal x,w,y,bottom,th; 845219b1045SBarry Smith ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 846219b1045SBarry Smith ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 8470298fd71SBarry Smith ierr = PetscDrawStringGetSize(draw,NULL,&th);CHKERRQ(ierr); 848219b1045SBarry Smith bottom = y - th; 849219b1045SBarry Smith for (i=levels-1; i>=0; i--) { 850b4375e8dSPeter Brune if (!mglevels[i]->smoothu || (mglevels[i]->smoothu == mglevels[i]->smoothd)) { 851219b1045SBarry Smith ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 852219b1045SBarry Smith ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 853219b1045SBarry Smith ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 854b4375e8dSPeter Brune } else { 855b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 856b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 857b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothd,viewer);CHKERRQ(ierr); 858b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 859b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 860b4375e8dSPeter Brune ierr = KSPView(mglevels[i]->smoothu,viewer);CHKERRQ(ierr); 861b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 862b4375e8dSPeter Brune } 8630298fd71SBarry Smith ierr = PetscDrawGetBoundingBox(draw,NULL,&bottom,NULL,NULL);CHKERRQ(ierr); 8641cd381d6SBarry Smith bottom -= th; 865219b1045SBarry Smith } 8665b0b0462SBarry Smith } 867f3fbd535SBarry Smith PetscFunctionReturn(0); 868f3fbd535SBarry Smith } 869f3fbd535SBarry Smith 870af0996ceSBarry Smith #include <petsc/private/kspimpl.h> 871cab2e9ccSBarry Smith 872f3fbd535SBarry Smith /* 873f3fbd535SBarry Smith Calls setup for the KSP on each level 874f3fbd535SBarry Smith */ 875f3fbd535SBarry Smith PetscErrorCode PCSetUp_MG(PC pc) 876f3fbd535SBarry Smith { 877f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 878f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 879f3fbd535SBarry Smith PetscErrorCode ierr; 8801a97d7b6SStefano Zampini PetscInt i,n; 88198e04984SBarry Smith PC cpc; 8823db492dfSBarry Smith PetscBool dump = PETSC_FALSE,opsset,use_amat,missinginterpolate = PETSC_FALSE; 883f3fbd535SBarry Smith Mat dA,dB; 884f3fbd535SBarry Smith Vec tvec; 885218a07d4SBarry Smith DM *dms; 8860a545947SLisandro Dalcin PetscViewer viewer = NULL; 88741b6fd38SMatthew G. Knepley PetscBool dAeqdB = PETSC_FALSE, needRestricts = PETSC_FALSE, doCR = PETSC_FALSE; 888f3fbd535SBarry Smith 889f3fbd535SBarry Smith PetscFunctionBegin; 8901a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels with PCMGSetLevels() before setting up"); 8911a97d7b6SStefano Zampini n = mglevels[0]->levels; 89201bc414fSDmitry Karpeev /* FIX: Move this to PCSetFromOptions_MG? */ 8933aeaf226SBarry Smith if (mg->usedmfornumberoflevels) { 8943aeaf226SBarry Smith PetscInt levels; 8953aeaf226SBarry Smith ierr = DMGetRefineLevel(pc->dm,&levels);CHKERRQ(ierr); 8963aeaf226SBarry Smith levels++; 8973aeaf226SBarry Smith if (levels > n) { /* the problem is now being solved on a finer grid */ 8980298fd71SBarry Smith ierr = PCMGSetLevels(pc,levels,NULL);CHKERRQ(ierr); 8993aeaf226SBarry Smith n = levels; 9003aeaf226SBarry Smith ierr = PCSetFromOptions(pc);CHKERRQ(ierr); /* it is bad to call this here, but otherwise will never be called for the new hierarchy */ 9013aeaf226SBarry Smith mglevels = mg->levels; 9023aeaf226SBarry Smith } 9033aeaf226SBarry Smith } 90498e04984SBarry Smith ierr = KSPGetPC(mglevels[0]->smoothd,&cpc);CHKERRQ(ierr); 9053aeaf226SBarry Smith 906f3fbd535SBarry Smith 907f3fbd535SBarry Smith /* If user did not provide fine grid operators OR operator was not updated since last global KSPSetOperators() */ 908f3fbd535SBarry Smith /* so use those from global PC */ 909f3fbd535SBarry Smith /* Is this what we always want? What if user wants to keep old one? */ 9100298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[n-1]->smoothd,NULL,&opsset);CHKERRQ(ierr); 91198e04984SBarry Smith if (opsset) { 91298e04984SBarry Smith Mat mmat; 91323ee1639SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,NULL,&mmat);CHKERRQ(ierr); 91498e04984SBarry Smith if (mmat == pc->pmat) opsset = PETSC_FALSE; 91598e04984SBarry Smith } 9164a5f07a5SMark F. Adams 91741b6fd38SMatthew G. Knepley /* Create CR solvers */ 91841b6fd38SMatthew G. Knepley ierr = PCMGGetAdaptCR(pc, &doCR);CHKERRQ(ierr); 91941b6fd38SMatthew G. Knepley if (doCR) { 92041b6fd38SMatthew G. Knepley const char *prefix; 92141b6fd38SMatthew G. Knepley 92241b6fd38SMatthew G. Knepley ierr = PCGetOptionsPrefix(pc, &prefix);CHKERRQ(ierr); 92341b6fd38SMatthew G. Knepley for (i = 1; i < n; ++i) { 92441b6fd38SMatthew G. Knepley PC ipc, cr; 92541b6fd38SMatthew G. Knepley char crprefix[128]; 92641b6fd38SMatthew G. Knepley 92741b6fd38SMatthew G. Knepley ierr = KSPCreate(PetscObjectComm((PetscObject) pc), &mglevels[i]->cr);CHKERRQ(ierr); 92841b6fd38SMatthew G. Knepley ierr = KSPSetErrorIfNotConverged(mglevels[i]->cr, PETSC_FALSE);CHKERRQ(ierr); 92941b6fd38SMatthew G. Knepley ierr = PetscObjectIncrementTabLevel((PetscObject) mglevels[i]->cr, (PetscObject) pc, n-i);CHKERRQ(ierr); 93041b6fd38SMatthew G. Knepley ierr = KSPSetOptionsPrefix(mglevels[i]->cr, prefix);CHKERRQ(ierr); 93141b6fd38SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) mglevels[i]->cr, PetscMGLevelId, mglevels[i]->level);CHKERRQ(ierr); 93241b6fd38SMatthew G. Knepley ierr = KSPSetType(mglevels[i]->cr, KSPCHEBYSHEV);CHKERRQ(ierr); 93341b6fd38SMatthew G. Knepley ierr = KSPSetConvergenceTest(mglevels[i]->cr, KSPConvergedSkip, NULL, NULL);CHKERRQ(ierr); 93441b6fd38SMatthew G. Knepley ierr = KSPSetNormType(mglevels[i]->cr, KSP_NORM_PRECONDITIONED);CHKERRQ(ierr); 93541b6fd38SMatthew G. Knepley ierr = KSPGetPC(mglevels[i]->cr, &ipc);CHKERRQ(ierr); 93641b6fd38SMatthew G. Knepley 93741b6fd38SMatthew G. Knepley ierr = PCSetType(ipc, PCCOMPOSITE);CHKERRQ(ierr); 93841b6fd38SMatthew G. Knepley ierr = PCCompositeSetType(ipc, PC_COMPOSITE_MULTIPLICATIVE);CHKERRQ(ierr); 93941b6fd38SMatthew G. Knepley ierr = PCCompositeAddPCType(ipc, PCSOR);CHKERRQ(ierr); 94041b6fd38SMatthew G. Knepley ierr = CreateCR_Private(pc, i, &cr);CHKERRQ(ierr); 94141b6fd38SMatthew G. Knepley ierr = PCCompositeAddPC(ipc, cr);CHKERRQ(ierr); 94241b6fd38SMatthew G. Knepley ierr = PCDestroy(&cr);CHKERRQ(ierr); 94341b6fd38SMatthew G. Knepley 94441b6fd38SMatthew G. Knepley ierr = KSPSetTolerances(mglevels[i]->cr, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, mg->default_smoothd);CHKERRQ(ierr); 94541b6fd38SMatthew G. Knepley ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr, PETSC_TRUE);CHKERRQ(ierr); 94641b6fd38SMatthew G. Knepley ierr = PetscSNPrintf(crprefix, 128, "mg_levels_%d_cr_", (int) i);CHKERRQ(ierr); 94741b6fd38SMatthew G. Knepley ierr = KSPAppendOptionsPrefix(mglevels[i]->cr, crprefix);CHKERRQ(ierr); 94841b6fd38SMatthew G. Knepley ierr = PetscLogObjectParent((PetscObject) pc, (PetscObject) mglevels[i]->cr);CHKERRQ(ierr); 94941b6fd38SMatthew G. Knepley } 95041b6fd38SMatthew G. Knepley } 95141b6fd38SMatthew G. Knepley 9524a5f07a5SMark F. Adams if (!opsset) { 95371b23a65SMark F. Adams ierr = PCGetUseAmat(pc,&use_amat);CHKERRQ(ierr); 9544a5f07a5SMark F. Adams if (use_amat) { 955f3fbd535SBarry 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); 95623ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->mat,pc->pmat);CHKERRQ(ierr); 95769858f1bSStefano Zampini } else { 9584a5f07a5SMark 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); 95923ee1639SBarry Smith ierr = KSPSetOperators(mglevels[n-1]->smoothd,pc->pmat,pc->pmat);CHKERRQ(ierr); 9604a5f07a5SMark F. Adams } 9614a5f07a5SMark F. Adams } 962f3fbd535SBarry Smith 9639c7ffb39SBarry Smith for (i=n-1; i>0; i--) { 9649c7ffb39SBarry Smith if (!(mglevels[i]->interpolate || mglevels[i]->restrct)) { 9659c7ffb39SBarry Smith missinginterpolate = PETSC_TRUE; 9669c7ffb39SBarry Smith continue; 9679c7ffb39SBarry Smith } 9689c7ffb39SBarry Smith } 9692134b1e4SBarry Smith 9702134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[n-1]->smoothd,&dA,&dB);CHKERRQ(ierr); 9712134b1e4SBarry Smith if (dA == dB) dAeqdB = PETSC_TRUE; 9722134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_NONE) || (((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_MAT)) && !dAeqdB)) { 9732134b1e4SBarry Smith needRestricts = PETSC_TRUE; /* user must compute either mat, pmat, or both so must restrict x to coarser levels */ 9742134b1e4SBarry Smith } 9752134b1e4SBarry Smith 9762134b1e4SBarry Smith 9779c7ffb39SBarry Smith /* 9789c7ffb39SBarry Smith Skipping if user has provided all interpolation/restriction needed (since DM might not be able to produce them (when coming from SNES/TS) 9799c7ffb39SBarry Smith Skipping for galerkin==2 (externally managed hierarchy such as ML and GAMG). Cleaner logic here would be great. Wrap ML/GAMG as DMs? 9809c7ffb39SBarry Smith */ 9812134b1e4SBarry Smith if (missinginterpolate && pc->dm && mg->galerkin != PC_MG_GALERKIN_EXTERNAL && !pc->setupcalled) { 9822d2b81a6SBarry Smith /* construct the interpolation from the DMs */ 9832e499ae9SBarry Smith Mat p; 98473250ac0SBarry Smith Vec rscale; 985785e854fSJed Brown ierr = PetscMalloc1(n,&dms);CHKERRQ(ierr); 986218a07d4SBarry Smith dms[n-1] = pc->dm; 987ef1267afSMatthew G. Knepley /* Separately create them so we do not get DMKSP interference between levels */ 988ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMCoarsen(dms[i+1],MPI_COMM_NULL,&dms[i]);CHKERRQ(ierr);} 98941fce8fdSHong Zhang /* 99041fce8fdSHong Zhang Force the mat type of coarse level operator to be AIJ because usually we want to use LU for coarse level. 99141fce8fdSHong Zhang Notice that it can be overwritten by -mat_type because KSPSetUp() reads command line options. 99241fce8fdSHong Zhang But it is safe to use -dm_mat_type. 99341fce8fdSHong Zhang 99441fce8fdSHong Zhang The mat type should not be hardcoded like this, we need to find a better way. 99541fce8fdSHong Zhang ierr = DMSetMatType(dms[0],MATAIJ);CHKERRQ(ierr); 99641fce8fdSHong Zhang */ 997218a07d4SBarry Smith for (i=n-2; i>-1; i--) { 998942e3340SBarry Smith DMKSP kdm; 999eab52d2bSLawrence Mitchell PetscBool dmhasrestrict, dmhasinject; 10003c0c59f3SBarry Smith ierr = KSPSetDM(mglevels[i]->smoothd,dms[i]);CHKERRQ(ierr); 10012134b1e4SBarry Smith if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothd,PETSC_FALSE);CHKERRQ(ierr);} 1002c27ee7a3SPatrick Farrell if (mglevels[i]->smoothd != mglevels[i]->smoothu) { 1003c27ee7a3SPatrick Farrell ierr = KSPSetDM(mglevels[i]->smoothu,dms[i]);CHKERRQ(ierr); 1004c27ee7a3SPatrick Farrell if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->smoothu,PETSC_FALSE);CHKERRQ(ierr);} 1005c27ee7a3SPatrick Farrell } 100641b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 100741b6fd38SMatthew G. Knepley ierr = KSPSetDM(mglevels[i]->cr,dms[i]);CHKERRQ(ierr); 100841b6fd38SMatthew G. Knepley if (!needRestricts) {ierr = KSPSetDMActive(mglevels[i]->cr,PETSC_FALSE);CHKERRQ(ierr);} 100941b6fd38SMatthew G. Knepley } 1010942e3340SBarry Smith ierr = DMGetDMKSPWrite(dms[i],&kdm);CHKERRQ(ierr); 1011d1a3e677SJed Brown /* Ugly hack so that the next KSPSetUp() will use the RHS that we set. A better fix is to change dmActive to take 1012d1a3e677SJed Brown * a bitwise OR of computing the matrix, RHS, and initial iterate. */ 10130298fd71SBarry Smith kdm->ops->computerhs = NULL; 10140298fd71SBarry Smith kdm->rhsctx = NULL; 101524c3aa18SBarry Smith if (!mglevels[i+1]->interpolate) { 1016e727c939SJed Brown ierr = DMCreateInterpolation(dms[i],dms[i+1],&p,&rscale);CHKERRQ(ierr); 10172d2b81a6SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,p);CHKERRQ(ierr); 101800ac8be1SBarry Smith if (rscale) {ierr = PCMGSetRScale(pc,i+1,rscale);CHKERRQ(ierr);} 101973250ac0SBarry Smith ierr = VecDestroy(&rscale);CHKERRQ(ierr); 10206bf464f9SBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 1021218a07d4SBarry Smith } 10223ad4599aSBarry Smith ierr = DMHasCreateRestriction(dms[i],&dmhasrestrict);CHKERRQ(ierr); 10233ad4599aSBarry Smith if (dmhasrestrict && !mglevels[i+1]->restrct){ 10243ad4599aSBarry Smith ierr = DMCreateRestriction(dms[i],dms[i+1],&p);CHKERRQ(ierr); 10253ad4599aSBarry Smith ierr = PCMGSetRestriction(pc,i+1,p);CHKERRQ(ierr); 10263ad4599aSBarry Smith ierr = MatDestroy(&p);CHKERRQ(ierr); 10273ad4599aSBarry Smith } 1028eab52d2bSLawrence Mitchell ierr = DMHasCreateInjection(dms[i],&dmhasinject);CHKERRQ(ierr); 1029eab52d2bSLawrence Mitchell if (dmhasinject && !mglevels[i+1]->inject){ 1030eab52d2bSLawrence Mitchell ierr = DMCreateInjection(dms[i],dms[i+1],&p);CHKERRQ(ierr); 1031eab52d2bSLawrence Mitchell ierr = PCMGSetInjection(pc,i+1,p);CHKERRQ(ierr); 1032eab52d2bSLawrence Mitchell ierr = MatDestroy(&p);CHKERRQ(ierr); 1033eab52d2bSLawrence Mitchell } 103424c3aa18SBarry Smith } 10352d2b81a6SBarry Smith 1036ef1267afSMatthew G. Knepley for (i=n-2; i>-1; i--) {ierr = DMDestroy(&dms[i]);CHKERRQ(ierr);} 10372d2b81a6SBarry Smith ierr = PetscFree(dms);CHKERRQ(ierr); 1038ce4cda84SJed Brown } 1039cab2e9ccSBarry Smith 1040ce4cda84SJed Brown if (pc->dm && !pc->setupcalled) { 10412134b1e4SBarry Smith /* finest smoother also gets DM but it is not active, independent of whether galerkin==PC_MG_GALERKIN_EXTERNAL */ 1042cab2e9ccSBarry Smith ierr = KSPSetDM(mglevels[n-1]->smoothd,pc->dm);CHKERRQ(ierr); 1043cab2e9ccSBarry Smith ierr = KSPSetDMActive(mglevels[n-1]->smoothd,PETSC_FALSE);CHKERRQ(ierr); 1044c27ee7a3SPatrick Farrell if (mglevels[n-1]->smoothd != mglevels[n-1]->smoothu) { 1045c27ee7a3SPatrick Farrell ierr = KSPSetDM(mglevels[n-1]->smoothu,pc->dm);CHKERRQ(ierr); 1046c27ee7a3SPatrick Farrell ierr = KSPSetDMActive(mglevels[n-1]->smoothu,PETSC_FALSE);CHKERRQ(ierr); 1047c27ee7a3SPatrick Farrell } 104841b6fd38SMatthew G. Knepley if (mglevels[n-1]->cr) { 104941b6fd38SMatthew G. Knepley ierr = KSPSetDM(mglevels[n-1]->cr,pc->dm);CHKERRQ(ierr); 105041b6fd38SMatthew G. Knepley ierr = KSPSetDMActive(mglevels[n-1]->cr,PETSC_FALSE);CHKERRQ(ierr); 105141b6fd38SMatthew G. Knepley } 1052218a07d4SBarry Smith } 1053218a07d4SBarry Smith 10542134b1e4SBarry Smith if (mg->galerkin < PC_MG_GALERKIN_NONE) { 10552134b1e4SBarry Smith Mat A,B; 10562134b1e4SBarry Smith PetscBool doA = PETSC_FALSE,doB = PETSC_FALSE; 10572134b1e4SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 10582134b1e4SBarry Smith 10592134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_PMAT) || (mg->galerkin == PC_MG_GALERKIN_BOTH)) doB = PETSC_TRUE; 10602134b1e4SBarry Smith if ((mg->galerkin == PC_MG_GALERKIN_MAT) || ((mg->galerkin == PC_MG_GALERKIN_BOTH) && (dA != dB))) doA = PETSC_TRUE; 10612134b1e4SBarry Smith if (pc->setupcalled) reuse = MAT_REUSE_MATRIX; 1062f3fbd535SBarry Smith for (i=n-2; i>-1; i--) { 1063b9367914SBarry 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"); 1064b9367914SBarry Smith if (!mglevels[i+1]->interpolate) { 1065b9367914SBarry Smith ierr = PCMGSetInterpolation(pc,i+1,mglevels[i+1]->restrct);CHKERRQ(ierr); 1066b9367914SBarry Smith } 1067b9367914SBarry Smith if (!mglevels[i+1]->restrct) { 1068b9367914SBarry Smith ierr = PCMGSetRestriction(pc,i+1,mglevels[i+1]->interpolate);CHKERRQ(ierr); 1069b9367914SBarry Smith } 10702134b1e4SBarry Smith if (reuse == MAT_REUSE_MATRIX) { 10712134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,&B);CHKERRQ(ierr); 10722134b1e4SBarry Smith } 10732134b1e4SBarry Smith if (doA) { 10742df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dA,mglevels[i+1]->interpolate,reuse,1.0,&A);CHKERRQ(ierr); 10752134b1e4SBarry Smith } 10762134b1e4SBarry Smith if (doB) { 10772df6c5c3SPatrick Farrell ierr = MatGalerkin(mglevels[i+1]->restrct,dB,mglevels[i+1]->interpolate,reuse,1.0,&B);CHKERRQ(ierr); 1078b9367914SBarry Smith } 10792134b1e4SBarry Smith /* the management of the PetscObjectReference() and PetscObjecDereference() below is rather delicate */ 10802134b1e4SBarry Smith if (!doA && dAeqdB) { 10812134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);} 10822134b1e4SBarry Smith A = B; 10832134b1e4SBarry Smith } else if (!doA && reuse == MAT_INITIAL_MATRIX) { 10842134b1e4SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&A,NULL);CHKERRQ(ierr); 10852134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1086b9367914SBarry Smith } 10872134b1e4SBarry Smith if (!doB && dAeqdB) { 10882134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) {ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);} 10892134b1e4SBarry Smith B = A; 10902134b1e4SBarry Smith } else if (!doB && reuse == MAT_INITIAL_MATRIX) { 109123ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,NULL,&B);CHKERRQ(ierr); 10922134b1e4SBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 10937d537d92SJed Brown } 10942134b1e4SBarry Smith if (reuse == MAT_INITIAL_MATRIX) { 10952134b1e4SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothd,A,B);CHKERRQ(ierr); 10962134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)A);CHKERRQ(ierr); 10972134b1e4SBarry Smith ierr = PetscObjectDereference((PetscObject)B);CHKERRQ(ierr); 10982134b1e4SBarry Smith } 10992134b1e4SBarry Smith dA = A; 1100f3fbd535SBarry Smith dB = B; 1101f3fbd535SBarry Smith } 1102f3fbd535SBarry Smith } 1103f3b08a26SMatthew G. Knepley 1104f3b08a26SMatthew G. Knepley 1105f3b08a26SMatthew G. Knepley /* Adapt interpolation matrices */ 1106f3b08a26SMatthew G. Knepley if (mg->adaptInterpolation) { 1107f3b08a26SMatthew G. Knepley mg->Nc = mg->Nc < 0 ? 6 : mg->Nc; /* Default to 6 modes */ 1108f3b08a26SMatthew G. Knepley for (i = 0; i < n; ++i) { 1109f3b08a26SMatthew G. Knepley ierr = PCMGComputeCoarseSpace_Internal(pc, i, mg->coarseSpaceType, mg->Nc, !i ? NULL : mglevels[i-1]->coarseSpace, &mglevels[i]->coarseSpace);CHKERRQ(ierr); 1110f3b08a26SMatthew 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);} 1111f3b08a26SMatthew G. Knepley } 1112f3b08a26SMatthew G. Knepley for (i = n-2; i > -1; --i) { 1113f3b08a26SMatthew G. Knepley ierr = PCMGRecomputeLevelOperators_Internal(pc, i);CHKERRQ(ierr); 1114f3b08a26SMatthew G. Knepley } 1115f3b08a26SMatthew G. Knepley } 1116f3b08a26SMatthew G. Knepley 11172134b1e4SBarry Smith if (needRestricts && pc->dm) { 1118caa4e7f2SJed Brown for (i=n-2; i>=0; i--) { 1119ccceb50cSJed Brown DM dmfine,dmcoarse; 1120caa4e7f2SJed Brown Mat Restrict,Inject; 1121caa4e7f2SJed Brown Vec rscale; 1122ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i+1]->smoothd,&dmfine);CHKERRQ(ierr); 1123ccceb50cSJed Brown ierr = KSPGetDM(mglevels[i]->smoothd,&dmcoarse);CHKERRQ(ierr); 1124caa4e7f2SJed Brown ierr = PCMGGetRestriction(pc,i+1,&Restrict);CHKERRQ(ierr); 1125caa4e7f2SJed Brown ierr = PCMGGetRScale(pc,i+1,&rscale);CHKERRQ(ierr); 1126eab52d2bSLawrence Mitchell ierr = PCMGGetInjection(pc,i+1,&Inject);CHKERRQ(ierr); 1127caa4e7f2SJed Brown ierr = DMRestrict(dmfine,Restrict,rscale,Inject,dmcoarse);CHKERRQ(ierr); 1128caa4e7f2SJed Brown } 1129caa4e7f2SJed Brown } 1130f3fbd535SBarry Smith 1131f3fbd535SBarry Smith if (!pc->setupcalled) { 1132f3fbd535SBarry Smith for (i=0; i<n; i++) { 1133f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothd);CHKERRQ(ierr); 1134f3fbd535SBarry Smith } 1135f3fbd535SBarry Smith for (i=1; i<n; i++) { 1136f3fbd535SBarry Smith if (mglevels[i]->smoothu && (mglevels[i]->smoothu != mglevels[i]->smoothd)) { 1137f3fbd535SBarry Smith ierr = KSPSetFromOptions(mglevels[i]->smoothu);CHKERRQ(ierr); 1138f3fbd535SBarry Smith } 113941b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 114041b6fd38SMatthew G. Knepley ierr = KSPSetFromOptions(mglevels[i]->cr);CHKERRQ(ierr); 114141b6fd38SMatthew G. Knepley } 1142f3fbd535SBarry Smith } 11433ad4599aSBarry Smith /* insure that if either interpolation or restriction is set the other other one is set */ 1144f3fbd535SBarry Smith for (i=1; i<n; i++) { 11453ad4599aSBarry Smith ierr = PCMGGetInterpolation(pc,i,NULL);CHKERRQ(ierr); 11463ad4599aSBarry Smith ierr = PCMGGetRestriction(pc,i,NULL);CHKERRQ(ierr); 1147f3fbd535SBarry Smith } 1148f3fbd535SBarry Smith for (i=0; i<n-1; i++) { 1149f3fbd535SBarry Smith if (!mglevels[i]->b) { 1150f3fbd535SBarry Smith Vec *vec; 11512a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[i]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 1152f3fbd535SBarry Smith ierr = PCMGSetRhs(pc,i,*vec);CHKERRQ(ierr); 11536bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 1154f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 1155f3fbd535SBarry Smith } 1156f3fbd535SBarry Smith if (!mglevels[i]->r && i) { 1157f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 1158f3fbd535SBarry Smith ierr = PCMGSetR(pc,i,tvec);CHKERRQ(ierr); 11596bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 1160f3fbd535SBarry Smith } 1161f3fbd535SBarry Smith if (!mglevels[i]->x) { 1162f3fbd535SBarry Smith ierr = VecDuplicate(mglevels[i]->b,&tvec);CHKERRQ(ierr); 1163f3fbd535SBarry Smith ierr = PCMGSetX(pc,i,tvec);CHKERRQ(ierr); 11646bf464f9SBarry Smith ierr = VecDestroy(&tvec);CHKERRQ(ierr); 1165f3fbd535SBarry Smith } 116641b6fd38SMatthew G. Knepley if (doCR) { 116741b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crx);CHKERRQ(ierr); 116841b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[i]->b,&mglevels[i]->crb);CHKERRQ(ierr); 116941b6fd38SMatthew G. Knepley ierr = VecZeroEntries(mglevels[i]->crb);CHKERRQ(ierr); 117041b6fd38SMatthew G. Knepley } 1171f3fbd535SBarry Smith } 1172f3fbd535SBarry Smith if (n != 1 && !mglevels[n-1]->r) { 1173f3fbd535SBarry Smith /* PCMGSetR() on the finest level if user did not supply it */ 1174f3fbd535SBarry Smith Vec *vec; 11752a7a6963SBarry Smith ierr = KSPCreateVecs(mglevels[n-1]->smoothd,1,&vec,0,NULL);CHKERRQ(ierr); 1176f3fbd535SBarry Smith ierr = PCMGSetR(pc,n-1,*vec);CHKERRQ(ierr); 11776bf464f9SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 1178f3fbd535SBarry Smith ierr = PetscFree(vec);CHKERRQ(ierr); 1179f3fbd535SBarry Smith } 118041b6fd38SMatthew G. Knepley if (doCR) { 118141b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crx);CHKERRQ(ierr); 118241b6fd38SMatthew G. Knepley ierr = VecDuplicate(mglevels[n-1]->r, &mglevels[n-1]->crb);CHKERRQ(ierr); 118341b6fd38SMatthew G. Knepley ierr = VecZeroEntries(mglevels[n-1]->crb);CHKERRQ(ierr); 118441b6fd38SMatthew G. Knepley } 1185f3fbd535SBarry Smith } 1186f3fbd535SBarry Smith 118798e04984SBarry Smith if (pc->dm) { 118898e04984SBarry Smith /* need to tell all the coarser levels to rebuild the matrix using the DM for that level */ 118998e04984SBarry Smith for (i=0; i<n-1; i++) { 119098e04984SBarry Smith if (mglevels[i]->smoothd->setupstage != KSP_SETUP_NEW) mglevels[i]->smoothd->setupstage = KSP_SETUP_NEWMATRIX; 119198e04984SBarry Smith } 119298e04984SBarry Smith } 1193f3fbd535SBarry Smith 1194f3fbd535SBarry Smith for (i=1; i<n; i++) { 11952a21e185SBarry Smith if (mglevels[i]->smoothu == mglevels[i]->smoothd || mg->am == PC_MG_FULL || mg->am == PC_MG_KASKADE || mg->cyclesperpcapply > 1){ 1196f3fbd535SBarry Smith /* if doing only down then initial guess is zero */ 1197f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothd,PETSC_TRUE);CHKERRQ(ierr); 1198f3fbd535SBarry Smith } 119941b6fd38SMatthew G. Knepley if (mglevels[i]->cr) {ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr);} 120063e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1201f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothd);CHKERRQ(ierr); 1202c0decd05SBarry Smith if (mglevels[i]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1203899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1204899639b0SHong Zhang } 120563e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1206d42688cbSBarry Smith if (!mglevels[i]->residual) { 1207d42688cbSBarry Smith Mat mat; 120813842ffbSBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 120954b2cd4bSJed Brown ierr = PCMGSetResidual(pc,i,PCMGResidualDefault,mat);CHKERRQ(ierr); 1210d42688cbSBarry Smith } 1211fcb023d4SJed Brown if (!mglevels[i]->residualtranspose) { 1212fcb023d4SJed Brown Mat mat; 1213fcb023d4SJed Brown ierr = KSPGetOperators(mglevels[i]->smoothd,&mat,NULL);CHKERRQ(ierr); 1214fcb023d4SJed Brown ierr = PCMGSetResidualTranspose(pc,i,PCMGResidualTransposeDefault,mat);CHKERRQ(ierr); 1215fcb023d4SJed Brown } 1216f3fbd535SBarry Smith } 1217f3fbd535SBarry Smith for (i=1; i<n; i++) { 1218f3fbd535SBarry Smith if (mglevels[i]->smoothu && mglevels[i]->smoothu != mglevels[i]->smoothd) { 1219f3fbd535SBarry Smith Mat downmat,downpmat; 1220f3fbd535SBarry Smith 1221f3fbd535SBarry Smith /* check if operators have been set for up, if not use down operators to set them */ 12220298fd71SBarry Smith ierr = KSPGetOperatorsSet(mglevels[i]->smoothu,&opsset,NULL);CHKERRQ(ierr); 1223f3fbd535SBarry Smith if (!opsset) { 122423ee1639SBarry Smith ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 122523ee1639SBarry Smith ierr = KSPSetOperators(mglevels[i]->smoothu,downmat,downpmat);CHKERRQ(ierr); 1226f3fbd535SBarry Smith } 1227f3fbd535SBarry Smith 1228f3fbd535SBarry Smith ierr = KSPSetInitialGuessNonzero(mglevels[i]->smoothu,PETSC_TRUE);CHKERRQ(ierr); 122963e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1230f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[i]->smoothu);CHKERRQ(ierr); 1231c0decd05SBarry Smith if (mglevels[i]->smoothu->reason == KSP_DIVERGED_PC_FAILED) { 1232899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1233899639b0SHong Zhang } 123463e6d426SJed Brown if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1235f3fbd535SBarry Smith } 123641b6fd38SMatthew G. Knepley if (mglevels[i]->cr) { 123741b6fd38SMatthew G. Knepley Mat downmat,downpmat; 123841b6fd38SMatthew G. Knepley 123941b6fd38SMatthew G. Knepley /* check if operators have been set for up, if not use down operators to set them */ 124041b6fd38SMatthew G. Knepley ierr = KSPGetOperatorsSet(mglevels[i]->cr,&opsset,NULL);CHKERRQ(ierr); 124141b6fd38SMatthew G. Knepley if (!opsset) { 124241b6fd38SMatthew G. Knepley ierr = KSPGetOperators(mglevels[i]->smoothd,&downmat,&downpmat);CHKERRQ(ierr); 124341b6fd38SMatthew G. Knepley ierr = KSPSetOperators(mglevels[i]->cr,downmat,downpmat);CHKERRQ(ierr); 124441b6fd38SMatthew G. Knepley } 124541b6fd38SMatthew G. Knepley 124641b6fd38SMatthew G. Knepley ierr = KSPSetInitialGuessNonzero(mglevels[i]->cr,PETSC_TRUE);CHKERRQ(ierr); 124741b6fd38SMatthew G. Knepley if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 124841b6fd38SMatthew G. Knepley ierr = KSPSetUp(mglevels[i]->cr);CHKERRQ(ierr); 124941b6fd38SMatthew G. Knepley if (mglevels[i]->cr->reason == KSP_DIVERGED_PC_FAILED) { 125041b6fd38SMatthew G. Knepley pc->failedreason = PC_SUBPC_ERROR; 125141b6fd38SMatthew G. Knepley } 125241b6fd38SMatthew G. Knepley if (mglevels[i]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 125341b6fd38SMatthew G. Knepley } 1254f3fbd535SBarry Smith } 1255f3fbd535SBarry Smith 125663e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventBegin(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1257f3fbd535SBarry Smith ierr = KSPSetUp(mglevels[0]->smoothd);CHKERRQ(ierr); 1258c0decd05SBarry Smith if (mglevels[0]->smoothd->reason == KSP_DIVERGED_PC_FAILED) { 1259899639b0SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 1260899639b0SHong Zhang } 126163e6d426SJed Brown if (mglevels[0]->eventsmoothsetup) {ierr = PetscLogEventEnd(mglevels[0]->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1262f3fbd535SBarry Smith 1263f3fbd535SBarry Smith /* 1264f3fbd535SBarry Smith Dump the interpolation/restriction matrices plus the 1265e3c5b3baSBarry Smith Jacobian/stiffness on each level. This allows MATLAB users to 1266f3fbd535SBarry Smith easily check if the Galerkin condition A_c = R A_f R^T is satisfied. 1267f3fbd535SBarry Smith 1268f3fbd535SBarry Smith Only support one or the other at the same time. 1269f3fbd535SBarry Smith */ 1270f3fbd535SBarry Smith #if defined(PETSC_USE_SOCKET_VIEWER) 1271c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_matlab",&dump,NULL);CHKERRQ(ierr); 1272ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_SOCKET_(PetscObjectComm((PetscObject)pc)); 1273f3fbd535SBarry Smith dump = PETSC_FALSE; 1274f3fbd535SBarry Smith #endif 1275c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)pc)->options,((PetscObject)pc)->prefix,"-pc_mg_dump_binary",&dump,NULL);CHKERRQ(ierr); 1276ce94432eSBarry Smith if (dump) viewer = PETSC_VIEWER_BINARY_(PetscObjectComm((PetscObject)pc)); 1277f3fbd535SBarry Smith 1278f3fbd535SBarry Smith if (viewer) { 1279f3fbd535SBarry Smith for (i=1; i<n; i++) { 1280f3fbd535SBarry Smith ierr = MatView(mglevels[i]->restrct,viewer);CHKERRQ(ierr); 1281f3fbd535SBarry Smith } 1282f3fbd535SBarry Smith for (i=0; i<n; i++) { 1283f3fbd535SBarry Smith ierr = KSPGetPC(mglevels[i]->smoothd,&pc);CHKERRQ(ierr); 1284f3fbd535SBarry Smith ierr = MatView(pc->mat,viewer);CHKERRQ(ierr); 1285f3fbd535SBarry Smith } 1286f3fbd535SBarry Smith } 1287f3fbd535SBarry Smith PetscFunctionReturn(0); 1288f3fbd535SBarry Smith } 1289f3fbd535SBarry Smith 1290f3fbd535SBarry Smith /* -------------------------------------------------------------------------------------*/ 1291f3fbd535SBarry Smith 129297d33e41SMatthew G. Knepley PetscErrorCode PCMGGetLevels_MG(PC pc, PetscInt *levels) 129397d33e41SMatthew G. Knepley { 129497d33e41SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 129597d33e41SMatthew G. Knepley 129697d33e41SMatthew G. Knepley PetscFunctionBegin; 129797d33e41SMatthew G. Knepley *levels = mg->nlevels; 129897d33e41SMatthew G. Knepley PetscFunctionReturn(0); 129997d33e41SMatthew G. Knepley } 130097d33e41SMatthew G. Knepley 13014b9ad928SBarry Smith /*@ 130297177400SBarry Smith PCMGGetLevels - Gets the number of levels to use with MG. 13034b9ad928SBarry Smith 13044b9ad928SBarry Smith Not Collective 13054b9ad928SBarry Smith 13064b9ad928SBarry Smith Input Parameter: 13074b9ad928SBarry Smith . pc - the preconditioner context 13084b9ad928SBarry Smith 13094b9ad928SBarry Smith Output parameter: 13104b9ad928SBarry Smith . levels - the number of levels 13114b9ad928SBarry Smith 13124b9ad928SBarry Smith Level: advanced 13134b9ad928SBarry Smith 131497177400SBarry Smith .seealso: PCMGSetLevels() 13154b9ad928SBarry Smith @*/ 13167087cfbeSBarry Smith PetscErrorCode PCMGGetLevels(PC pc,PetscInt *levels) 13174b9ad928SBarry Smith { 1318e952c529SMatthew G. Knepley PetscErrorCode ierr; 13194b9ad928SBarry Smith 13204b9ad928SBarry Smith PetscFunctionBegin; 13210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13224482741eSBarry Smith PetscValidIntPointer(levels,2); 132397d33e41SMatthew G. Knepley *levels = 0; 132437b5128cSMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGGetLevels_C",(PC,PetscInt*),(pc,levels));CHKERRQ(ierr); 13254b9ad928SBarry Smith PetscFunctionReturn(0); 13264b9ad928SBarry Smith } 13274b9ad928SBarry Smith 13284b9ad928SBarry Smith /*@ 132997177400SBarry Smith PCMGSetType - Determines the form of multigrid to use: 13304b9ad928SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 13314b9ad928SBarry Smith 1332ad4df100SBarry Smith Logically Collective on PC 13334b9ad928SBarry Smith 13344b9ad928SBarry Smith Input Parameters: 13354b9ad928SBarry Smith + pc - the preconditioner context 13369dcbbd2bSBarry Smith - form - multigrid form, one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE, 13379dcbbd2bSBarry Smith PC_MG_FULL, PC_MG_KASKADE 13384b9ad928SBarry Smith 13394b9ad928SBarry Smith Options Database Key: 13404b9ad928SBarry Smith . -pc_mg_type <form> - Sets <form>, one of multiplicative, 13414b9ad928SBarry Smith additive, full, kaskade 13424b9ad928SBarry Smith 13434b9ad928SBarry Smith Level: advanced 13444b9ad928SBarry Smith 134597177400SBarry Smith .seealso: PCMGSetLevels() 13464b9ad928SBarry Smith @*/ 13477087cfbeSBarry Smith PetscErrorCode PCMGSetType(PC pc,PCMGType form) 13484b9ad928SBarry Smith { 1349f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 13504b9ad928SBarry Smith 13514b9ad928SBarry Smith PetscFunctionBegin; 13520700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1353c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,form,2); 135431567311SBarry Smith mg->am = form; 13559dcbbd2bSBarry Smith if (form == PC_MG_MULTIPLICATIVE) pc->ops->applyrichardson = PCApplyRichardson_MG; 1356c60c7ad4SBarry Smith else pc->ops->applyrichardson = NULL; 1357c60c7ad4SBarry Smith PetscFunctionReturn(0); 1358c60c7ad4SBarry Smith } 1359c60c7ad4SBarry Smith 1360c60c7ad4SBarry Smith /*@ 1361c60c7ad4SBarry Smith PCMGGetType - Determines the form of multigrid to use: 1362c60c7ad4SBarry Smith multiplicative, additive, full, or the Kaskade algorithm. 1363c60c7ad4SBarry Smith 1364c60c7ad4SBarry Smith Logically Collective on PC 1365c60c7ad4SBarry Smith 1366c60c7ad4SBarry Smith Input Parameter: 1367c60c7ad4SBarry Smith . pc - the preconditioner context 1368c60c7ad4SBarry Smith 1369c60c7ad4SBarry Smith Output Parameter: 1370c60c7ad4SBarry Smith . type - one of PC_MG_MULTIPLICATIVE, PC_MG_ADDITIVE,PC_MG_FULL, PC_MG_KASKADE 1371c60c7ad4SBarry Smith 1372c60c7ad4SBarry Smith 1373c60c7ad4SBarry Smith Level: advanced 1374c60c7ad4SBarry Smith 1375c60c7ad4SBarry Smith .seealso: PCMGSetLevels() 1376c60c7ad4SBarry Smith @*/ 1377c60c7ad4SBarry Smith PetscErrorCode PCMGGetType(PC pc,PCMGType *type) 1378c60c7ad4SBarry Smith { 1379c60c7ad4SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1380c60c7ad4SBarry Smith 1381c60c7ad4SBarry Smith PetscFunctionBegin; 1382c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1383c60c7ad4SBarry Smith *type = mg->am; 13844b9ad928SBarry Smith PetscFunctionReturn(0); 13854b9ad928SBarry Smith } 13864b9ad928SBarry Smith 13874b9ad928SBarry Smith /*@ 13880d353602SBarry Smith PCMGSetCycleType - Sets the type cycles to use. Use PCMGSetCycleTypeOnLevel() for more 13894b9ad928SBarry Smith complicated cycling. 13904b9ad928SBarry Smith 1391ad4df100SBarry Smith Logically Collective on PC 13924b9ad928SBarry Smith 13934b9ad928SBarry Smith Input Parameters: 1394c2be2410SBarry Smith + pc - the multigrid context 1395c1cbb1deSBarry Smith - n - either PC_MG_CYCLE_V or PC_MG_CYCLE_W 13964b9ad928SBarry Smith 13974b9ad928SBarry Smith Options Database Key: 1398c1cbb1deSBarry Smith . -pc_mg_cycle_type <v,w> - provide the cycle desired 13994b9ad928SBarry Smith 14004b9ad928SBarry Smith Level: advanced 14014b9ad928SBarry Smith 14020d353602SBarry Smith .seealso: PCMGSetCycleTypeOnLevel() 14034b9ad928SBarry Smith @*/ 14047087cfbeSBarry Smith PetscErrorCode PCMGSetCycleType(PC pc,PCMGCycleType n) 14054b9ad928SBarry Smith { 1406f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1407f3fbd535SBarry Smith PC_MG_Levels **mglevels = mg->levels; 140879416396SBarry Smith PetscInt i,levels; 14094b9ad928SBarry Smith 14104b9ad928SBarry Smith PetscFunctionBegin; 14110700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 141269fbec6eSBarry Smith PetscValidLogicalCollectiveEnum(pc,n,2); 14131a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1414f3fbd535SBarry Smith levels = mglevels[0]->levels; 14152fa5cd67SKarl Rupp for (i=0; i<levels; i++) mglevels[i]->cycles = n; 14164b9ad928SBarry Smith PetscFunctionReturn(0); 14174b9ad928SBarry Smith } 14184b9ad928SBarry Smith 14198cc2d5dfSBarry Smith /*@ 14208cc2d5dfSBarry Smith PCMGMultiplicativeSetCycles - Sets the number of cycles to use for each preconditioner step 14218cc2d5dfSBarry Smith of multigrid when PCMGType of PC_MG_MULTIPLICATIVE is used 14228cc2d5dfSBarry Smith 1423ad4df100SBarry Smith Logically Collective on PC 14248cc2d5dfSBarry Smith 14258cc2d5dfSBarry Smith Input Parameters: 14268cc2d5dfSBarry Smith + pc - the multigrid context 14278cc2d5dfSBarry Smith - n - number of cycles (default is 1) 14288cc2d5dfSBarry Smith 14298cc2d5dfSBarry Smith Options Database Key: 1430e1bc860dSBarry Smith . -pc_mg_multiplicative_cycles n 14318cc2d5dfSBarry Smith 14328cc2d5dfSBarry Smith Level: advanced 14338cc2d5dfSBarry Smith 143495452b02SPatrick Sanan Notes: 143595452b02SPatrick Sanan This is not associated with setting a v or w cycle, that is set with PCMGSetCycleType() 14368cc2d5dfSBarry Smith 14378cc2d5dfSBarry Smith .seealso: PCMGSetCycleTypeOnLevel(), PCMGSetCycleType() 14388cc2d5dfSBarry Smith @*/ 14397087cfbeSBarry Smith PetscErrorCode PCMGMultiplicativeSetCycles(PC pc,PetscInt n) 14408cc2d5dfSBarry Smith { 1441f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 14428cc2d5dfSBarry Smith 14438cc2d5dfSBarry Smith PetscFunctionBegin; 14440700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1445c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 14462a21e185SBarry Smith mg->cyclesperpcapply = n; 14478cc2d5dfSBarry Smith PetscFunctionReturn(0); 14488cc2d5dfSBarry Smith } 14498cc2d5dfSBarry Smith 14502134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin_MG(PC pc,PCMGGalerkinType use) 1451fb15c04eSBarry Smith { 1452fb15c04eSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1453fb15c04eSBarry Smith 1454fb15c04eSBarry Smith PetscFunctionBegin; 14552134b1e4SBarry Smith mg->galerkin = use; 1456fb15c04eSBarry Smith PetscFunctionReturn(0); 1457fb15c04eSBarry Smith } 1458fb15c04eSBarry Smith 1459c2be2410SBarry Smith /*@ 146097177400SBarry Smith PCMGSetGalerkin - Causes the coarser grid matrices to be computed from the 146182b87d87SMatthew G. Knepley finest grid via the Galerkin process: A_i-1 = r_i * A_i * p_i 1462c2be2410SBarry Smith 1463ad4df100SBarry Smith Logically Collective on PC 1464c2be2410SBarry Smith 1465c2be2410SBarry Smith Input Parameters: 1466c91913d3SJed Brown + pc - the multigrid context 14672134b1e4SBarry Smith - use - one of PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, or PC_MG_GALERKIN_NONE 1468c2be2410SBarry Smith 1469c2be2410SBarry Smith Options Database Key: 14702134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> 1471c2be2410SBarry Smith 1472c2be2410SBarry Smith Level: intermediate 1473c2be2410SBarry Smith 147495452b02SPatrick Sanan Notes: 147595452b02SPatrick Sanan Some codes that use PCMG such as PCGAMG use Galerkin internally while constructing the hierarchy and thus do not 14761c1aac46SBarry Smith use the PCMG construction of the coarser grids. 14771c1aac46SBarry Smith 14782134b1e4SBarry Smith .seealso: PCMGGetGalerkin(), PCMGGalerkinType 14793fc8bf9cSBarry Smith 1480c2be2410SBarry Smith @*/ 14812134b1e4SBarry Smith PetscErrorCode PCMGSetGalerkin(PC pc,PCMGGalerkinType use) 1482c2be2410SBarry Smith { 1483fb15c04eSBarry Smith PetscErrorCode ierr; 1484c2be2410SBarry Smith 1485c2be2410SBarry Smith PetscFunctionBegin; 14860700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 14872134b1e4SBarry Smith ierr = PetscTryMethod(pc,"PCMGSetGalerkin_C",(PC,PCMGGalerkinType),(pc,use));CHKERRQ(ierr); 1488c2be2410SBarry Smith PetscFunctionReturn(0); 1489c2be2410SBarry Smith } 1490c2be2410SBarry Smith 14913fc8bf9cSBarry Smith /*@ 14923fc8bf9cSBarry Smith PCMGGetGalerkin - Checks if Galerkin multigrid is being used, i.e. 149382b87d87SMatthew G. Knepley A_i-1 = r_i * A_i * p_i 14943fc8bf9cSBarry Smith 14953fc8bf9cSBarry Smith Not Collective 14963fc8bf9cSBarry Smith 14973fc8bf9cSBarry Smith Input Parameter: 14983fc8bf9cSBarry Smith . pc - the multigrid context 14993fc8bf9cSBarry Smith 15003fc8bf9cSBarry Smith Output Parameter: 15012134b1e4SBarry 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 15023fc8bf9cSBarry Smith 15033fc8bf9cSBarry Smith Level: intermediate 15043fc8bf9cSBarry Smith 15052134b1e4SBarry Smith .seealso: PCMGSetGalerkin(), PCMGGalerkinType 15063fc8bf9cSBarry Smith 15073fc8bf9cSBarry Smith @*/ 15082134b1e4SBarry Smith PetscErrorCode PCMGGetGalerkin(PC pc,PCMGGalerkinType *galerkin) 15093fc8bf9cSBarry Smith { 1510f3fbd535SBarry Smith PC_MG *mg = (PC_MG*)pc->data; 15113fc8bf9cSBarry Smith 15123fc8bf9cSBarry Smith PetscFunctionBegin; 15130700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 15142134b1e4SBarry Smith *galerkin = mg->galerkin; 15153fc8bf9cSBarry Smith PetscFunctionReturn(0); 15163fc8bf9cSBarry Smith } 15173fc8bf9cSBarry Smith 1518f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation_MG(PC pc, PetscBool adapt) 1519f3b08a26SMatthew G. Knepley { 1520f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1521f3b08a26SMatthew G. Knepley 1522f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1523f3b08a26SMatthew G. Knepley mg->adaptInterpolation = adapt; 1524f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1525f3b08a26SMatthew G. Knepley } 1526f3b08a26SMatthew G. Knepley 1527f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation_MG(PC pc, PetscBool *adapt) 1528f3b08a26SMatthew G. Knepley { 1529f3b08a26SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 1530f3b08a26SMatthew G. Knepley 1531f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1532f3b08a26SMatthew G. Knepley *adapt = mg->adaptInterpolation; 1533f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1534f3b08a26SMatthew G. Knepley } 1535f3b08a26SMatthew G. Knepley 153641b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR_MG(PC pc, PetscBool cr) 153741b6fd38SMatthew G. Knepley { 153841b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 153941b6fd38SMatthew G. Knepley 154041b6fd38SMatthew G. Knepley PetscFunctionBegin; 154141b6fd38SMatthew G. Knepley mg->compatibleRelaxation = cr; 154241b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 154341b6fd38SMatthew G. Knepley } 154441b6fd38SMatthew G. Knepley 154541b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR_MG(PC pc, PetscBool *cr) 154641b6fd38SMatthew G. Knepley { 154741b6fd38SMatthew G. Knepley PC_MG *mg = (PC_MG *) pc->data; 154841b6fd38SMatthew G. Knepley 154941b6fd38SMatthew G. Knepley PetscFunctionBegin; 155041b6fd38SMatthew G. Knepley *cr = mg->compatibleRelaxation; 155141b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 155241b6fd38SMatthew G. Knepley } 155341b6fd38SMatthew G. Knepley 1554f3b08a26SMatthew G. Knepley /*@ 1555f3b08a26SMatthew 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. 1556f3b08a26SMatthew G. Knepley 1557f3b08a26SMatthew G. Knepley Logically Collective on PC 1558f3b08a26SMatthew G. Knepley 1559f3b08a26SMatthew G. Knepley Input Parameters: 1560f3b08a26SMatthew G. Knepley + pc - the multigrid context 1561f3b08a26SMatthew G. Knepley - adapt - flag for adaptation of the interpolator 1562f3b08a26SMatthew G. Knepley 1563f3b08a26SMatthew G. Knepley Options Database Keys: 1564f3b08a26SMatthew G. Knepley + -pc_mg_adapt_interp - Turn on adaptation 1565f3b08a26SMatthew G. Knepley . -pc_mg_adapt_interp_n <int> - The number of modes to use, should be divisible by dimension 1566f3b08a26SMatthew G. Knepley - -pc_mg_adapt_interp_coarse_space <type> - The type of coarse space: polynomial, harmonic, eigenvector, generalized_eigenvector 1567f3b08a26SMatthew G. Knepley 1568f3b08a26SMatthew G. Knepley Level: intermediate 1569f3b08a26SMatthew G. Knepley 1570f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1571f3b08a26SMatthew G. Knepley .seealso: PCMGGetAdaptInterpolation(), PCMGSetGalerkin() 1572f3b08a26SMatthew G. Knepley @*/ 1573f3b08a26SMatthew G. Knepley PetscErrorCode PCMGSetAdaptInterpolation(PC pc, PetscBool adapt) 1574f3b08a26SMatthew G. Knepley { 1575f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1576f3b08a26SMatthew G. Knepley 1577f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1578f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1579f3b08a26SMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetAdaptInterpolation_C",(PC,PetscBool),(pc,adapt));CHKERRQ(ierr); 1580f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1581f3b08a26SMatthew G. Knepley } 1582f3b08a26SMatthew G. Knepley 1583f3b08a26SMatthew G. Knepley /*@ 1584f3b08a26SMatthew 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. 1585f3b08a26SMatthew G. Knepley 1586f3b08a26SMatthew G. Knepley Not collective 1587f3b08a26SMatthew G. Knepley 1588f3b08a26SMatthew G. Knepley Input Parameter: 1589f3b08a26SMatthew G. Knepley . pc - the multigrid context 1590f3b08a26SMatthew G. Knepley 1591f3b08a26SMatthew G. Knepley Output Parameter: 1592f3b08a26SMatthew G. Knepley . adapt - flag for adaptation of the interpolator 1593f3b08a26SMatthew G. Knepley 1594f3b08a26SMatthew G. Knepley Level: intermediate 1595f3b08a26SMatthew G. Knepley 1596f3b08a26SMatthew G. Knepley .keywords: MG, set, Galerkin 1597f3b08a26SMatthew G. Knepley .seealso: PCMGSetAdaptInterpolation(), PCMGSetGalerkin() 1598f3b08a26SMatthew G. Knepley @*/ 1599f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetAdaptInterpolation(PC pc, PetscBool *adapt) 1600f3b08a26SMatthew G. Knepley { 1601f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1602f3b08a26SMatthew G. Knepley 1603f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1604f3b08a26SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1605f3b08a26SMatthew G. Knepley PetscValidPointer(adapt, 2); 1606f3b08a26SMatthew G. Knepley ierr = PetscUseMethod(pc,"PCMGGetAdaptInterpolation_C",(PC,PetscBool*),(pc,adapt));CHKERRQ(ierr); 1607f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1608f3b08a26SMatthew G. Knepley } 1609f3b08a26SMatthew G. Knepley 16104b9ad928SBarry Smith /*@ 161141b6fd38SMatthew G. Knepley PCMGSetAdaptCR - Monitor the coarse space quality using an auxiliary solve with compatible relaxation. 161241b6fd38SMatthew G. Knepley 161341b6fd38SMatthew G. Knepley Logically Collective on PC 161441b6fd38SMatthew G. Knepley 161541b6fd38SMatthew G. Knepley Input Parameters: 161641b6fd38SMatthew G. Knepley + pc - the multigrid context 161741b6fd38SMatthew G. Knepley - cr - flag for compatible relaxation 161841b6fd38SMatthew G. Knepley 161941b6fd38SMatthew G. Knepley Options Database Keys: 162041b6fd38SMatthew G. Knepley . -pc_mg_adapt_cr - Turn on compatible relaxation 162141b6fd38SMatthew G. Knepley 162241b6fd38SMatthew G. Knepley Level: intermediate 162341b6fd38SMatthew G. Knepley 162441b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin 162541b6fd38SMatthew G. Knepley .seealso: PCMGGetAdaptCR(), PCMGSetAdaptInterpolation(), PCMGSetGalerkin() 162641b6fd38SMatthew G. Knepley @*/ 162741b6fd38SMatthew G. Knepley PetscErrorCode PCMGSetAdaptCR(PC pc, PetscBool cr) 162841b6fd38SMatthew G. Knepley { 162941b6fd38SMatthew G. Knepley PetscErrorCode ierr; 163041b6fd38SMatthew G. Knepley 163141b6fd38SMatthew G. Knepley PetscFunctionBegin; 163241b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 163341b6fd38SMatthew G. Knepley ierr = PetscTryMethod(pc,"PCMGSetAdaptCR_C",(PC,PetscBool),(pc,cr));CHKERRQ(ierr); 163441b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 163541b6fd38SMatthew G. Knepley } 163641b6fd38SMatthew G. Knepley 163741b6fd38SMatthew G. Knepley /*@ 163841b6fd38SMatthew G. Knepley PCMGGetAdaptCR - Get the flag to monitor coarse space quality using an auxiliary solve with compatible relaxation. 163941b6fd38SMatthew G. Knepley 164041b6fd38SMatthew G. Knepley Not collective 164141b6fd38SMatthew G. Knepley 164241b6fd38SMatthew G. Knepley Input Parameter: 164341b6fd38SMatthew G. Knepley . pc - the multigrid context 164441b6fd38SMatthew G. Knepley 164541b6fd38SMatthew G. Knepley Output Parameter: 164641b6fd38SMatthew G. Knepley . cr - flag for compatible relaxaion 164741b6fd38SMatthew G. Knepley 164841b6fd38SMatthew G. Knepley Level: intermediate 164941b6fd38SMatthew G. Knepley 165041b6fd38SMatthew G. Knepley .keywords: MG, set, Galerkin 165141b6fd38SMatthew G. Knepley .seealso: PCMGSetAdaptCR(), PCMGGetAdaptInterpolation(), PCMGSetGalerkin() 165241b6fd38SMatthew G. Knepley @*/ 165341b6fd38SMatthew G. Knepley PetscErrorCode PCMGGetAdaptCR(PC pc, PetscBool *cr) 165441b6fd38SMatthew G. Knepley { 165541b6fd38SMatthew G. Knepley PetscErrorCode ierr; 165641b6fd38SMatthew G. Knepley 165741b6fd38SMatthew G. Knepley PetscFunctionBegin; 165841b6fd38SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 165941b6fd38SMatthew G. Knepley PetscValidPointer(cr, 2); 166041b6fd38SMatthew G. Knepley ierr = PetscUseMethod(pc,"PCMGGetAdaptCR_C",(PC,PetscBool*),(pc,cr));CHKERRQ(ierr); 166141b6fd38SMatthew G. Knepley PetscFunctionReturn(0); 166241b6fd38SMatthew G. Knepley } 166341b6fd38SMatthew G. Knepley 166441b6fd38SMatthew G. Knepley /*@ 166506792cafSBarry Smith PCMGSetNumberSmooth - Sets the number of pre and post-smoothing steps to use 1666710315b6SLawrence Mitchell on all levels. Use PCMGDistinctSmoothUp() to create separate up and down smoothers if you want different numbers of 1667710315b6SLawrence Mitchell pre- and post-smoothing steps. 166806792cafSBarry Smith 166906792cafSBarry Smith Logically Collective on PC 167006792cafSBarry Smith 167106792cafSBarry Smith Input Parameters: 167206792cafSBarry Smith + mg - the multigrid context 167306792cafSBarry Smith - n - the number of smoothing steps 167406792cafSBarry Smith 167506792cafSBarry Smith Options Database Key: 1676a2b725a8SWilliam Gropp . -mg_levels_ksp_max_it <n> - Sets number of pre and post-smoothing steps 167706792cafSBarry Smith 167806792cafSBarry Smith Level: advanced 167906792cafSBarry Smith 168095452b02SPatrick Sanan Notes: 168195452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 168206792cafSBarry Smith there is no separate smooth up on the coarsest grid. 168306792cafSBarry Smith 1684710315b6SLawrence Mitchell .seealso: PCMGSetDistinctSmoothUp() 168506792cafSBarry Smith @*/ 168606792cafSBarry Smith PetscErrorCode PCMGSetNumberSmooth(PC pc,PetscInt n) 168706792cafSBarry Smith { 168806792cafSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 168906792cafSBarry Smith PC_MG_Levels **mglevels = mg->levels; 169006792cafSBarry Smith PetscErrorCode ierr; 169106792cafSBarry Smith PetscInt i,levels; 169206792cafSBarry Smith 169306792cafSBarry Smith PetscFunctionBegin; 169406792cafSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 169506792cafSBarry Smith PetscValidLogicalCollectiveInt(pc,n,2); 16961a97d7b6SStefano Zampini if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 169706792cafSBarry Smith levels = mglevels[0]->levels; 169806792cafSBarry Smith 169906792cafSBarry Smith for (i=1; i<levels; i++) { 170006792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothu,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 170106792cafSBarry Smith ierr = KSPSetTolerances(mglevels[i]->smoothd,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,n);CHKERRQ(ierr); 170206792cafSBarry Smith mg->default_smoothu = n; 170306792cafSBarry Smith mg->default_smoothd = n; 170406792cafSBarry Smith } 170506792cafSBarry Smith PetscFunctionReturn(0); 170606792cafSBarry Smith } 170706792cafSBarry Smith 1708f442ab6aSBarry Smith /*@ 17098e5aa403SBarry Smith PCMGSetDistinctSmoothUp - sets the up (post) smoother to be a separate KSP from the down (pre) smoother on all levels 1710710315b6SLawrence Mitchell and adds the suffix _up to the options name 1711f442ab6aSBarry Smith 1712f442ab6aSBarry Smith Logically Collective on PC 1713f442ab6aSBarry Smith 1714f442ab6aSBarry Smith Input Parameters: 1715f442ab6aSBarry Smith . pc - the preconditioner context 1716f442ab6aSBarry Smith 1717f442ab6aSBarry Smith Options Database Key: 1718f442ab6aSBarry Smith . -pc_mg_distinct_smoothup 1719f442ab6aSBarry Smith 1720f442ab6aSBarry Smith Level: advanced 1721f442ab6aSBarry Smith 172295452b02SPatrick Sanan Notes: 172395452b02SPatrick Sanan this does not set a value on the coarsest grid, since we assume that 1724f442ab6aSBarry Smith there is no separate smooth up on the coarsest grid. 1725f442ab6aSBarry Smith 1726710315b6SLawrence Mitchell .seealso: PCMGSetNumberSmooth() 1727f442ab6aSBarry Smith @*/ 1728f442ab6aSBarry Smith PetscErrorCode PCMGSetDistinctSmoothUp(PC pc) 1729f442ab6aSBarry Smith { 1730f442ab6aSBarry Smith PC_MG *mg = (PC_MG*)pc->data; 1731f442ab6aSBarry Smith PC_MG_Levels **mglevels = mg->levels; 1732f442ab6aSBarry Smith PetscErrorCode ierr; 1733f442ab6aSBarry Smith PetscInt i,levels; 1734f442ab6aSBarry Smith KSP subksp; 1735f442ab6aSBarry Smith 1736f442ab6aSBarry Smith PetscFunctionBegin; 1737f442ab6aSBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1738f442ab6aSBarry Smith if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Must set MG levels with PCMGSetLevels() before calling"); 1739f442ab6aSBarry Smith levels = mglevels[0]->levels; 1740f442ab6aSBarry Smith 1741f442ab6aSBarry Smith for (i=1; i<levels; i++) { 1742710315b6SLawrence Mitchell const char *prefix = NULL; 1743f442ab6aSBarry Smith /* make sure smoother up and down are different */ 1744f442ab6aSBarry Smith ierr = PCMGGetSmootherUp(pc,i,&subksp);CHKERRQ(ierr); 1745710315b6SLawrence Mitchell ierr = KSPGetOptionsPrefix(mglevels[i]->smoothd,&prefix);CHKERRQ(ierr); 1746710315b6SLawrence Mitchell ierr = KSPSetOptionsPrefix(subksp,prefix);CHKERRQ(ierr); 1747f442ab6aSBarry Smith ierr = KSPAppendOptionsPrefix(subksp,"up_");CHKERRQ(ierr); 1748f442ab6aSBarry Smith } 1749f442ab6aSBarry Smith PetscFunctionReturn(0); 1750f442ab6aSBarry Smith } 1751f442ab6aSBarry Smith 175207a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 175307a4832bSFande Kong PetscErrorCode PCGetInterpolations_MG(PC pc,PetscInt *num_levels,Mat *interpolations[]) 175407a4832bSFande Kong { 175507a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 175607a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 175707a4832bSFande Kong Mat *mat; 175807a4832bSFande Kong PetscInt l; 175907a4832bSFande Kong PetscErrorCode ierr; 176007a4832bSFande Kong 176107a4832bSFande Kong PetscFunctionBegin; 176207a4832bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 176307a4832bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 176407a4832bSFande Kong for (l=1; l< mg->nlevels; l++) { 176507a4832bSFande Kong mat[l-1] = mglevels[l]->interpolate; 176607a4832bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l-1]);CHKERRQ(ierr); 176707a4832bSFande Kong } 176807a4832bSFande Kong *num_levels = mg->nlevels; 176907a4832bSFande Kong *interpolations = mat; 177007a4832bSFande Kong PetscFunctionReturn(0); 177107a4832bSFande Kong } 177207a4832bSFande Kong 177307a4832bSFande Kong /* No new matrices are created, and the coarse operator matrices are the references to the original ones */ 177407a4832bSFande Kong PetscErrorCode PCGetCoarseOperators_MG(PC pc,PetscInt *num_levels,Mat *coarseOperators[]) 177507a4832bSFande Kong { 177607a4832bSFande Kong PC_MG *mg = (PC_MG*)pc->data; 177707a4832bSFande Kong PC_MG_Levels **mglevels = mg->levels; 177807a4832bSFande Kong PetscInt l; 177907a4832bSFande Kong Mat *mat; 178007a4832bSFande Kong PetscErrorCode ierr; 178107a4832bSFande Kong 178207a4832bSFande Kong PetscFunctionBegin; 178307a4832bSFande Kong if (!mglevels) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set MG levels before calling"); 178407a4832bSFande Kong ierr = PetscMalloc1(mg->nlevels,&mat);CHKERRQ(ierr); 178507a4832bSFande Kong for (l=0; l<mg->nlevels-1; l++) { 178607a4832bSFande Kong ierr = KSPGetOperators(mglevels[l]->smoothd,NULL,&(mat[l]));CHKERRQ(ierr); 178707a4832bSFande Kong ierr = PetscObjectReference((PetscObject)mat[l]);CHKERRQ(ierr); 178807a4832bSFande Kong } 178907a4832bSFande Kong *num_levels = mg->nlevels; 179007a4832bSFande Kong *coarseOperators = mat; 179107a4832bSFande Kong PetscFunctionReturn(0); 179207a4832bSFande Kong } 179307a4832bSFande Kong 1794f3b08a26SMatthew G. Knepley /*@C 1795f3b08a26SMatthew G. Knepley PCMGRegisterCoarseSpaceConstructor - Adds a method to the PCMG package for coarse space construction. 1796f3b08a26SMatthew G. Knepley 1797f3b08a26SMatthew G. Knepley Not collective 1798f3b08a26SMatthew G. Knepley 1799f3b08a26SMatthew G. Knepley Input Parameters: 1800f3b08a26SMatthew G. Knepley + name - name of the constructor 1801f3b08a26SMatthew G. Knepley - function - constructor routine 1802f3b08a26SMatthew G. Knepley 1803f3b08a26SMatthew G. Knepley Notes: 1804f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1805f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp) 1806f3b08a26SMatthew G. Knepley $ pc - The PC object 1807f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1808f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1809f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1810f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1811f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1812f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1813f3b08a26SMatthew G. Knepley 1814f3b08a26SMatthew G. Knepley Level: advanced 1815f3b08a26SMatthew G. Knepley 1816f3b08a26SMatthew G. Knepley .seealso: PCMGGetCoarseSpaceConstructor(), PCRegister() 1817f3b08a26SMatthew G. Knepley @*/ 1818f3b08a26SMatthew G. Knepley PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char name[], PetscErrorCode (*function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)) 1819f3b08a26SMatthew G. Knepley { 1820f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1821f3b08a26SMatthew G. Knepley 1822f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1823f3b08a26SMatthew G. Knepley ierr = PCInitializePackage();CHKERRQ(ierr); 1824f3b08a26SMatthew G. Knepley ierr = PetscFunctionListAdd(&PCMGCoarseList,name,function);CHKERRQ(ierr); 1825f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1826f3b08a26SMatthew G. Knepley } 1827f3b08a26SMatthew G. Knepley 1828f3b08a26SMatthew G. Knepley /*@C 1829f3b08a26SMatthew G. Knepley PCMGGetCoarseSpaceConstructor - Returns the given coarse space construction method. 1830f3b08a26SMatthew G. Knepley 1831f3b08a26SMatthew G. Knepley Not collective 1832f3b08a26SMatthew G. Knepley 1833f3b08a26SMatthew G. Knepley Input Parameter: 1834f3b08a26SMatthew G. Knepley . name - name of the constructor 1835f3b08a26SMatthew G. Knepley 1836f3b08a26SMatthew G. Knepley Output Parameter: 1837f3b08a26SMatthew G. Knepley . function - constructor routine 1838f3b08a26SMatthew G. Knepley 1839f3b08a26SMatthew G. Knepley Notes: 1840f3b08a26SMatthew G. Knepley Calling sequence for the routine: 1841f3b08a26SMatthew G. Knepley $ my_csp(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, const Vec initGuess[], Vec **coarseSp) 1842f3b08a26SMatthew G. Knepley $ pc - The PC object 1843f3b08a26SMatthew G. Knepley $ l - The multigrid level, 0 is the coarse level 1844f3b08a26SMatthew G. Knepley $ dm - The DM for this level 1845f3b08a26SMatthew G. Knepley $ smooth - The level smoother 1846f3b08a26SMatthew G. Knepley $ Nc - The size of the coarse space 1847f3b08a26SMatthew G. Knepley $ initGuess - Basis for an initial guess for the space 1848f3b08a26SMatthew G. Knepley $ coarseSp - A basis for the computed coarse space 1849f3b08a26SMatthew G. Knepley 1850f3b08a26SMatthew G. Knepley Level: advanced 1851f3b08a26SMatthew G. Knepley 1852f3b08a26SMatthew G. Knepley .seealso: PCMGRegisterCoarseSpaceConstructor(), PCRegister() 1853f3b08a26SMatthew G. Knepley @*/ 1854f3b08a26SMatthew G. Knepley PetscErrorCode PCMGGetCoarseSpaceConstructor(const char name[], PetscErrorCode (**function)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)) 1855f3b08a26SMatthew G. Knepley { 1856f3b08a26SMatthew G. Knepley PetscErrorCode ierr; 1857f3b08a26SMatthew G. Knepley 1858f3b08a26SMatthew G. Knepley PetscFunctionBegin; 1859f3b08a26SMatthew G. Knepley ierr = PetscFunctionListFind(PCMGCoarseList,name,function);CHKERRQ(ierr); 1860f3b08a26SMatthew G. Knepley PetscFunctionReturn(0); 1861f3b08a26SMatthew G. Knepley } 1862f3b08a26SMatthew G. Knepley 18634b9ad928SBarry Smith /* ----------------------------------------------------------------------------------------*/ 18644b9ad928SBarry Smith 18653b09bd56SBarry Smith /*MC 1866ccb205f8SBarry Smith PCMG - Use multigrid preconditioning. This preconditioner requires you provide additional 18673b09bd56SBarry Smith information about the coarser grid matrices and restriction/interpolation operators. 18683b09bd56SBarry Smith 18693b09bd56SBarry Smith Options Database Keys: 18703b09bd56SBarry Smith + -pc_mg_levels <nlevels> - number of levels including finest 1871391689abSStefano Zampini . -pc_mg_cycle_type <v,w> - provide the cycle desired 18728c1c2452SJed Brown . -pc_mg_type <additive,multiplicative,full,kaskade> - multiplicative is the default 18733b09bd56SBarry Smith . -pc_mg_log - log information about time spent on each level of the solver 1874710315b6SLawrence Mitchell . -pc_mg_distinct_smoothup - configure up (after interpolation) and down (before restriction) smoothers separately (with different options prefixes) 18752134b1e4SBarry Smith . -pc_mg_galerkin <both,pmat,mat,none> - use Galerkin process to compute coarser operators, i.e. Acoarse = R A R' 18768cf6037eSBarry Smith . -pc_mg_multiplicative_cycles - number of cycles to use as the preconditioner (defaults to 1) 18778cf6037eSBarry Smith . -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 1878e3c5b3baSBarry Smith to the Socket viewer for reading from MATLAB. 18798cf6037eSBarry Smith - -pc_mg_dump_binary - dumps the matrices for each level and the restriction/interpolation matrices 18808cf6037eSBarry Smith to the binary output file called binaryoutput 18813b09bd56SBarry Smith 188295452b02SPatrick Sanan Notes: 18830b3c753eSRichard 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 18843b09bd56SBarry Smith 18858cf6037eSBarry Smith When run with a single level the smoother options are used on that level NOT the coarse grid solver options 18868cf6037eSBarry Smith 188723067569SBarry Smith When run with KSPRICHARDSON the convergence test changes slightly if monitor is turned on. The iteration count may change slightly. This 188823067569SBarry Smith is because without monitoring the residual norm is computed WITHIN each multigrid cycle on the finest level after the pre-smoothing 188923067569SBarry Smith (because the residual has just been computed for the multigrid algorithm and is hence available for free) while with monitoring the 189023067569SBarry Smith residual is computed at the end of each cycle. 189123067569SBarry Smith 18923b09bd56SBarry Smith Level: intermediate 18933b09bd56SBarry Smith 18948cf6037eSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, PCEXOTIC, PCGAMG, PCML, PCHYPRE 1895710315b6SLawrence Mitchell PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), PCMGSetCycleType(), 1896710315b6SLawrence Mitchell PCMGSetDistinctSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 189797177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 18980d353602SBarry Smith PCMGSetCycleTypeOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 18993b09bd56SBarry Smith M*/ 19003b09bd56SBarry Smith 19018cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_MG(PC pc) 19024b9ad928SBarry Smith { 1903f3fbd535SBarry Smith PC_MG *mg; 1904f3fbd535SBarry Smith PetscErrorCode ierr; 1905f3fbd535SBarry Smith 19064b9ad928SBarry Smith PetscFunctionBegin; 1907b00a9115SJed Brown ierr = PetscNewLog(pc,&mg);CHKERRQ(ierr); 1908f3fbd535SBarry Smith pc->data = (void*)mg; 1909f3fbd535SBarry Smith mg->nlevels = -1; 191010eca3edSLisandro Dalcin mg->am = PC_MG_MULTIPLICATIVE; 19112134b1e4SBarry Smith mg->galerkin = PC_MG_GALERKIN_NONE; 1912f3b08a26SMatthew G. Knepley mg->adaptInterpolation = PETSC_FALSE; 1913f3b08a26SMatthew G. Knepley mg->Nc = -1; 1914f3b08a26SMatthew G. Knepley mg->eigenvalue = -1; 1915f3fbd535SBarry Smith 191637a44384SMark Adams pc->useAmat = PETSC_TRUE; 191737a44384SMark Adams 19184b9ad928SBarry Smith pc->ops->apply = PCApply_MG; 1919fcb023d4SJed Brown pc->ops->applytranspose = PCApplyTranspose_MG; 1920*30b0564aSStefano Zampini pc->ops->matapply = PCMatApply_MG; 19214b9ad928SBarry Smith pc->ops->setup = PCSetUp_MG; 1922a06653b4SBarry Smith pc->ops->reset = PCReset_MG; 19234b9ad928SBarry Smith pc->ops->destroy = PCDestroy_MG; 19244b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_MG; 19254b9ad928SBarry Smith pc->ops->view = PCView_MG; 1926fb15c04eSBarry Smith 1927f3b08a26SMatthew G. Knepley ierr = PetscObjectComposedDataRegister(&mg->eigenvalue);CHKERRQ(ierr); 1928fb15c04eSBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetGalerkin_C",PCMGSetGalerkin_MG);CHKERRQ(ierr); 192997d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetLevels_C",PCMGGetLevels_MG);CHKERRQ(ierr); 193097d33e41SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetLevels_C",PCMGSetLevels_MG);CHKERRQ(ierr); 1931fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetInterpolations_C",PCGetInterpolations_MG);CHKERRQ(ierr); 1932fd2dd295SFande Kong ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGetCoarseOperators_C",PCGetCoarseOperators_MG);CHKERRQ(ierr); 1933f3b08a26SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptInterpolation_C",PCMGSetAdaptInterpolation_MG);CHKERRQ(ierr); 1934f3b08a26SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptInterpolation_C",PCMGGetAdaptInterpolation_MG);CHKERRQ(ierr); 193541b6fd38SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGSetAdaptCR_C",PCMGSetAdaptCR_MG);CHKERRQ(ierr); 193641b6fd38SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCMGGetAdaptCR_C",PCMGGetAdaptCR_MG);CHKERRQ(ierr); 19374b9ad928SBarry Smith PetscFunctionReturn(0); 19384b9ad928SBarry Smith } 1939