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