static PetscErrorCode PCApply_MG(PC pc,Vec b,Vec x) { PC_MG *mg = (PC_MG*)pc->data; PC_MG_Levels **mglevels = mg->levels; PetscErrorCode ierr; PetscInt levels = mglevels[0]->levels,i; PetscFunctionBegin; if (mg->stageApply) {ierr = PetscLogStagePush(mg->stageApply);CHKERRQ(ierr);} /* When the DM is supplying the matrix then it will not exist until here */ for (i=0; i<levels; i++) { if (!mglevels[i]->A) { ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); } } mglevels[levels-1]->b = b; mglevels[levels-1]->x = x; if (mg->am == PC_MG_MULTIPLICATIVE) { ierr = VecSet(x,0.0);CHKERRQ(ierr); for (i=0; i<mg->cyclesperpcapply; i++) { ierr = PCMGMCycle_Private(pc,mglevels+levels-1,NULL);CHKERRQ(ierr); } } else if (mg->am == PC_MG_ADDITIVE) { ierr = PCMGACycle_Private(pc,mglevels);CHKERRQ(ierr); } else if (mg->am == PC_MG_KASKADE) { ierr = PCMGKCycle_Private(pc,mglevels);CHKERRQ(ierr); } else { ierr = PCMGFCycle_Private(pc,mglevels);CHKERRQ(ierr); } if (mg->stageApply) {ierr = PetscLogStagePop();CHKERRQ(ierr);} PetscFunctionReturn(0); }
PetscErrorCode PCMGMCycle_Private(PC pc,PC_MG_Levels **mglevelsin,PCRichardsonConvergedReason *reason) { PC_MG *mg = (PC_MG*)pc->data; PC_MG_Levels *mgc,*mglevels = *mglevelsin; PetscErrorCode ierr; PetscInt cycles = (mglevels->level == 1) ? 1 : (PetscInt) mglevels->cycles; PC subpc; PCFailedReason pcreason; PetscFunctionBegin; if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} ierr = KSPSolve(mglevels->smoothd,mglevels->b,mglevels->x);CHKERRQ(ierr); /* pre-smooth */ ierr = KSPGetPC(mglevels->smoothd,&subpc);CHKERRQ(ierr); ierr = PCGetSetUpFailedReason(subpc,&pcreason);CHKERRQ(ierr); if (pcreason) { pc->failedreason = PC_SUBPC_ERROR; } if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} if (mglevels->level) { /* not the coarsest grid */ if (mglevels->eventresidual) {ierr = PetscLogEventBegin(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} ierr = (*mglevels->residual)(mglevels->A,mglevels->b,mglevels->x,mglevels->r);CHKERRQ(ierr); if (mglevels->eventresidual) {ierr = PetscLogEventEnd(mglevels->eventresidual,0,0,0,0);CHKERRQ(ierr);} /* if on finest level and have convergence criteria set */ if (mglevels->level == mglevels->levels-1 && mg->ttol && reason) { PetscReal rnorm; ierr = VecNorm(mglevels->r,NORM_2,&rnorm);CHKERRQ(ierr); if (rnorm <= mg->ttol) { if (rnorm < mg->abstol) { *reason = PCRICHARDSON_CONVERGED_ATOL; ierr = PetscInfo2(pc,"Linear solver has converged. Residual norm %g is less than absolute tolerance %g\n",(double)rnorm,(double)mg->abstol);CHKERRQ(ierr); } else { *reason = PCRICHARDSON_CONVERGED_RTOL; 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); } PetscFunctionReturn(0); } } mgc = *(mglevelsin - 1); if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} ierr = MatRestrict(mglevels->restrct,mglevels->r,mgc->b);CHKERRQ(ierr); if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} ierr = VecSet(mgc->x,0.0);CHKERRQ(ierr); while (cycles--) { ierr = PCMGMCycle_Private(pc,mglevelsin-1,reason);CHKERRQ(ierr); } if (mglevels->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} ierr = MatInterpolateAdd(mglevels->interpolate,mgc->x,mglevels->x,mglevels->x);CHKERRQ(ierr); if (mglevels->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} if (mglevels->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} ierr = KSPSolve(mglevels->smoothu,mglevels->b,mglevels->x);CHKERRQ(ierr); /* post smooth */ if (mglevels->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} } PetscFunctionReturn(0); }
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) { PC_MG *mg = (PC_MG*)pc->data; PC_MG_Levels **mglevels = mg->levels; PetscErrorCode ierr; PetscInt levels = mglevels[0]->levels,i; PetscFunctionBegin; /* When the DM is supplying the matrix then it will not exist until here */ for (i=0; i<levels; i++) { if (!mglevels[i]->A) { ierr = KSPGetOperators(mglevels[i]->smoothu,&mglevels[i]->A,NULL);CHKERRQ(ierr); ierr = PetscObjectReference((PetscObject)mglevels[i]->A);CHKERRQ(ierr); } } mglevels[levels-1]->b = b; mglevels[levels-1]->x = x; mg->rtol = rtol; mg->abstol = abstol; mg->dtol = dtol; if (rtol) { /* compute initial residual norm for relative convergence test */ PetscReal rnorm; if (zeroguess) { ierr = VecNorm(b,NORM_2,&rnorm);CHKERRQ(ierr); } else { ierr = (*mglevels[levels-1]->residual)(mglevels[levels-1]->A,b,x,w);CHKERRQ(ierr); ierr = VecNorm(w,NORM_2,&rnorm);CHKERRQ(ierr); } mg->ttol = PetscMax(rtol*rnorm,abstol); } else if (abstol) mg->ttol = abstol; else mg->ttol = 0.0; /* since smoother is applied to full system, not just residual we need to make sure that smoothers don't stop prematurely due to small residual */ for (i=1; i<levels; i++) { ierr = KSPSetTolerances(mglevels[i]->smoothu,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); if (mglevels[i]->smoothu != mglevels[i]->smoothd) { ierr = KSPSetTolerances(mglevels[i]->smoothd,0,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); } } *reason = (PCRichardsonConvergedReason)0; for (i=0; i<its; i++) { ierr = PCMGMCycle_Private(pc,mglevels+levels-1,reason);CHKERRQ(ierr); if (*reason) break; } if (!*reason) *reason = PCRICHARDSON_CONVERGED_ITS; *outits = i; PetscFunctionReturn(0); }
PetscErrorCode PCMGFCycle_Private(PC pc,PC_MG_Levels **mglevels) { PetscErrorCode ierr; PetscInt i,l = mglevels[0]->levels; PetscFunctionBegin; /* restrict the RHS through all levels to coarsest. */ for (i=l-1; i>0; i--){ if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr); if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} } /* work our way up through the levels */ ierr = VecSet(mglevels[0]->x,0.0);CHKERRQ(ierr); for (i=0; i<l-1; i++) { ierr = PCMGMCycle_Private(pc,&mglevels[i],PETSC_NULL);CHKERRQ(ierr); if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} ierr = MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);CHKERRQ(ierr); if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} } ierr = PCMGMCycle_Private(pc,&mglevels[l-1],PETSC_NULL);CHKERRQ(ierr); PetscFunctionReturn(0); }