xref: /petsc/src/ksp/pc/impls/mg/fmg.c (revision 906ed7cc33fecbafab01746eec64dcdcc8a4842f)
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 PCMGMCycle_Private(PC_MG **,PetscTruth*);
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "PCMGFCycle_Private"
10 PetscErrorCode PCMGFCycle_Private(PC_MG **mg)
11 {
12   PetscErrorCode ierr;
13   PetscInt       i,l = mg[0]->levels;
14 
15   PetscFunctionBegin;
16   /* restrict the RHS through all levels to coarsest. */
17   for (i=l-1; i>0; i--){
18     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
19   }
20 
21   /* work our way up through the levels */
22   ierr = VecSet(mg[0]->x,0.0);CHKERRQ(ierr);
23   for (i=0; i<l-1; i++) {
24     ierr = PCMGMCycle_Private(&mg[i],PETSC_NULL);CHKERRQ(ierr);
25     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
26   }
27   ierr = PCMGMCycle_Private(&mg[l-1],PETSC_NULL);CHKERRQ(ierr);
28   PetscFunctionReturn(0);
29 }
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "PCMGKCycle_Private"
33 PetscErrorCode PCMGKCycle_Private(PC_MG **mg)
34 {
35   PetscErrorCode ierr;
36   PetscInt       i,l = mg[0]->levels;
37 
38   PetscFunctionBegin;
39   /* restrict the RHS through all levels to coarsest. */
40   for (i=l-1; i>0; i--){
41     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
42   }
43 
44   /* work our way up through the levels */
45   ierr = VecSet(mg[0]->x,0.0);CHKERRQ(ierr);
46   for (i=0; i<l-1; i++) {
47     if (mg[i]->eventsolve) {ierr = PetscLogEventBegin(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
48     ierr = KSPSolve(mg[i]->smoothd,mg[i]->b,mg[i]->x);CHKERRQ(ierr);
49     if (mg[i]->eventsolve) {ierr = PetscLogEventEnd(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
50     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
51   }
52   if (mg[l-1]->eventsolve) {ierr = PetscLogEventBegin(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
53   ierr = KSPSolve(mg[l-1]->smoothd,mg[l-1]->b,mg[l-1]->x);CHKERRQ(ierr);
54   if (mg[l-1]->eventsolve) {ierr = PetscLogEventEnd(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
55 
56   PetscFunctionReturn(0);
57 }
58 
59 
60