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