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