xref: /petsc/src/ksp/pc/impls/mg/fmg.c (revision 702e6c57188f6ca39e32e26fd8fcae05ffbfed3e)
1 /*$Id: fmg.c,v 1.26 2001/08/21 21:03:20 bsmith Exp $*/
2 /*
3      Full multigrid using either additive or multiplicative V or W cycle
4 */
5 #include "src/ksp/pc/impls/mg/mgimpl.h"
6 
7 EXTERN int MGMCycle_Private(MG *,PetscTruth*);
8 
9 /*
10        MGFCycle_Private - Given an MG structure created with MGCreate() runs
11                full multigrid.
12 
13     Iput Parameters:
14 .   mg - structure created with MGCreate().
15 
16     Note: This may not be what others call full multigrid. What we
17           do is restrict the rhs to all levels, then starting
18           on the coarsest level work our way up generating
19           initial guess for the next level. This provides an
20           improved preconditioner but not a great improvement.
21 */
22 #undef __FUNCT__
23 #define __FUNCT__ "MGFCycle_Private"
24 int MGFCycle_Private(MG *mg)
25 {
26   int    i,l = mg[0]->levels,ierr;
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 int MGKCycle_Private(MG *mg)
57 {
58   int    i,l = mg[0]->levels,ierr;
59   PetscScalar zero = 0.0;
60 
61   PetscFunctionBegin;
62   /* restrict the RHS through all levels to coarsest. */
63   for (i=l-1; i>0; i--){
64     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
65   }
66 
67   /* work our way up through the levels */
68   ierr = VecSet(&zero,mg[0]->x);CHKERRQ(ierr);
69   for (i=0; i<l-1; i++) {
70     if (mg[i]->eventsolve) {ierr = PetscLogEventBegin(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
71     ierr = KSPSetRhs(mg[i]->smoothd,mg[i]->b);CHKERRQ(ierr);
72     ierr = KSPSetSolution(mg[i]->smoothd,mg[i]->x);CHKERRQ(ierr);
73     ierr = KSPSolve(mg[i]->smoothd);CHKERRQ(ierr);
74     if (mg[i]->eventsolve) {ierr = PetscLogEventEnd(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
75     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
76   }
77   if (mg[l-1]->eventsolve) {ierr = PetscLogEventBegin(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
78   ierr = KSPSetRhs(mg[l-1]->smoothd,mg[l-1]->b);CHKERRQ(ierr);
79   ierr = KSPSetSolution(mg[l-1]->smoothd,mg[l-1]->x);CHKERRQ(ierr);
80   ierr = KSPSolve(mg[l-1]->smoothd);CHKERRQ(ierr);
81   if (mg[l-1]->eventsolve) {ierr = PetscLogEventEnd(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
82 
83   PetscFunctionReturn(0);
84 }
85 
86 
87