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