Actual source code: smg.c
petsc-3.6.1 2015-08-06
2: /*
3: Additive Multigrid V Cycle routine
4: */
5: #include <petsc/private/pcmgimpl.h>
9: PetscErrorCode PCMGACycle_Private(PC pc,PC_MG_Levels **mglevels)
10: {
12: PetscInt i,l = mglevels[0]->levels;
15: /* compute RHS on each level */
16: for (i=l-1; i>0; i--) {
17: if (mglevels[i]->eventinterprestrict) {PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);}
18: MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);
19: if (mglevels[i]->eventinterprestrict) {PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);}
20: }
21: /* solve separately on each level */
22: for (i=0; i<l; i++) {
23: VecSet(mglevels[i]->x,0.0);
24: if (mglevels[i]->eventsmoothsolve) {PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0);}
25: KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);
26: if (mglevels[i]->eventsmoothsolve) {PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0);}
27: }
28: for (i=1; i<l; i++) {
29: if (mglevels[i]->eventinterprestrict) {PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);}
30: MatInterpolateAdd(mglevels[i]->interpolate,mglevels[i-1]->x,mglevels[i]->x,mglevels[i]->x);
31: if (mglevels[i]->eventinterprestrict) {PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);}
32: }
33: return(0);
34: }