xref: /petsc/src/ksp/pc/impls/mg/fmg.c (revision ecd4d7eeb1f43d54a7849690752177b3c01f4bdc)
1 /*
2      Full multigrid using either additive or multiplicative V or W cycle
3 */
4 #include "src/ksp/pc/impls/mg/mgimpl.h"
5 
6 EXTERN PetscErrorCode MGMCycle_Private(MG *,PetscTruth*);
7 
8 /*
9        MGFCycle_Private - Given an MG structure created with MGCreate() runs
10                full multigrid.
11 
12     Iput Parameters:
13 .   mg - structure created with MGCreate().
14 
15     Note: This may not be what others call full multigrid. What we
16           do is restrict the rhs to all levels, then starting
17           on the coarsest level work our way up generating
18           initial guess for the next level. This provides an
19           improved preconditioner but not a great improvement.
20 */
21 #undef __FUNCT__
22 #define __FUNCT__ "MGFCycle_Private"
23 PetscErrorCode MGFCycle_Private(MG *mg)
24 {
25   PetscErrorCode ierr;
26   int    i,l = mg[0]->levels;
27   PetscScalar zero = 0.0;
28 
29   PetscFunctionBegin;
30   /* restrict the RHS through all levels to coarsest. */
31   for (i=l-1; i>0; i--){
32     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
33   }
34 
35   /* work our way up through the levels */
36   ierr = VecSet(&zero,mg[0]->x);CHKERRQ(ierr);
37   for (i=0; i<l-1; i++) {
38     ierr = MGMCycle_Private(&mg[i],PETSC_NULL);CHKERRQ(ierr);
39     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
40   }
41   ierr = MGMCycle_Private(&mg[l-1],PETSC_NULL);CHKERRQ(ierr);
42   PetscFunctionReturn(0);
43 }
44 
45 /*
46        MGKCycle_Private - Given an MG structure created with MGCreate() runs
47                full Kascade MG solve.
48 
49     Iput Parameters:
50 .   mg - structure created with MGCreate().
51 
52     Note: This may not be what others call Kascadic MG.
53 */
54 #undef __FUNCT__
55 #define __FUNCT__ "MGKCycle_Private"
56 PetscErrorCode MGKCycle_Private(MG *mg)
57 {
58   PetscErrorCode ierr;
59   int    i,l = mg[0]->levels;
60   PetscScalar zero = 0.0;
61 
62   PetscFunctionBegin;
63   /* restrict the RHS through all levels to coarsest. */
64   for (i=l-1; i>0; i--){
65     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
66   }
67 
68   /* work our way up through the levels */
69   ierr = VecSet(&zero,mg[0]->x);CHKERRQ(ierr);
70   for (i=0; i<l-1; i++) {
71     if (mg[i]->eventsolve) {ierr = PetscLogEventBegin(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
72     ierr = KSPSolve(mg[i]->smoothd,mg[i]->b,mg[i]->x);CHKERRQ(ierr);
73     if (mg[i]->eventsolve) {ierr = PetscLogEventEnd(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
74     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
75   }
76   if (mg[l-1]->eventsolve) {ierr = PetscLogEventBegin(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
77   ierr = KSPSolve(mg[l-1]->smoothd,mg[l-1]->b,mg[l-1]->x);CHKERRQ(ierr);
78   if (mg[l-1]->eventsolve) {ierr = PetscLogEventEnd(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
79 
80   PetscFunctionReturn(0);
81 }
82 
83 
84