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