xref: /petsc/src/ksp/pc/impls/mg/smg.c (revision 5e71baeff2f3138f93cd4f5927dfd596eb8325cc)
1 
2 /*
3      Additive Multigrid V Cycle routine
4 */
5 #include <petsc/private/pcmgimpl.h>
6 
7 PetscErrorCode PCMGACycle_Private(PC pc,PC_MG_Levels **mglevels)
8 {
9   PetscErrorCode ierr;
10   PetscInt       i,l = mglevels[0]->levels;
11 
12   PetscFunctionBegin;
13   /* compute RHS on each level */
14   for (i=l-1; i>0; i--) {
15     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
16     ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr);
17     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
18   }
19   /* solve separately on each level */
20   for (i=0; i<l; i++) {
21     ierr = VecSet(mglevels[i]->x,0.0);CHKERRQ(ierr);
22     if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
23     ierr = KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);CHKERRQ(ierr);
24     if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
25   }
26   for (i=1; i<l; i++) {
27     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
28     ierr = MatInterpolateAdd(mglevels[i]->interpolate,mglevels[i-1]->x,mglevels[i]->x,mglevels[i]->x);CHKERRQ(ierr);
29     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
30   }
31   PetscFunctionReturn(0);
32 }
33