1: /*
2: Full multigrid using either additive or multiplicative V or W cycle
3: */
4: #include <petsc/private/pcmgimpl.h>
6: extern PetscErrorCode PCMGMCycle_Private(PC,PC_MG_Levels**,PCRichardsonConvergedReason*);
10: PetscErrorCode PCMGFCycle_Private(PC pc,PC_MG_Levels **mglevels) 11: {
13: PetscInt i,l = mglevels[0]->levels;
16: /* restrict the RHS through all levels to coarsest. */
17: for (i=l-1; i>0; i--) {
18: if (mglevels[i]->eventinterprestrict) {PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);}
19: MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);
20: if (mglevels[i]->eventinterprestrict) {PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);}
21: }
23: /* work our way up through the levels */
24: VecSet(mglevels[0]->x,0.0);
25: for (i=0; i<l-1; i++) {
26: PCMGMCycle_Private(pc,&mglevels[i],NULL);
27: if (mglevels[i+1]->eventinterprestrict) {PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);}
28: MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);
29: if (mglevels[i+1]->eventinterprestrict) {PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);}
30: }
31: PCMGMCycle_Private(pc,&mglevels[l-1],NULL);
32: return(0);
33: }
37: PetscErrorCode PCMGKCycle_Private(PC pc,PC_MG_Levels **mglevels) 38: {
40: PetscInt i,l = mglevels[0]->levels;
43: /* restrict the RHS through all levels to coarsest. */
44: for (i=l-1; i>0; i--) {
45: if (mglevels[i]->eventinterprestrict) {PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);}
46: MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);
47: if (mglevels[i]->eventinterprestrict) {PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);}
48: }
50: /* work our way up through the levels */
51: VecSet(mglevels[0]->x,0.0);
52: for (i=0; i<l-1; i++) {
53: if (mglevels[i]->eventsmoothsolve) {PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0);}
54: KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);
55: if (mglevels[i]->eventsmoothsolve) {PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0);}
56: if (mglevels[i+1]->eventinterprestrict) {PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);}
57: MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);
58: if (mglevels[i+1]->eventinterprestrict) {PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);}
59: }
60: if (mglevels[l-1]->eventsmoothsolve) {PetscLogEventBegin(mglevels[l-1]->eventsmoothsolve,0,0,0,0);}
61: KSPSolve(mglevels[l-1]->smoothd,mglevels[l-1]->b,mglevels[l-1]->x);
62: if (mglevels[l-1]->eventsmoothsolve) {PetscLogEventEnd(mglevels[l-1]->eventsmoothsolve,0,0,0,0);}
63: return(0);
64: }