xref: /petsc/src/ksp/pc/impls/mg/smg.c (revision 794163968b27cf3436b159913eecff082c02cff0)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith      Additive Multigrid V Cycle routine
34b9ad928SBarry Smith */
44b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h"
54b9ad928SBarry Smith 
64b9ad928SBarry Smith /*
74b9ad928SBarry Smith        MGACycle_Private - Given an MG structure created with MGCreate() runs
84b9ad928SBarry Smith                   one cycle down through the levels and back up. Applys
94b9ad928SBarry Smith                   the smoothers in an additive manner.
104b9ad928SBarry Smith 
114b9ad928SBarry Smith     Iput Parameters:
124b9ad928SBarry Smith .   mg - structure created with  MGCreate().
134b9ad928SBarry Smith 
144b9ad928SBarry Smith */
154b9ad928SBarry Smith #undef __FUNCT__
164b9ad928SBarry Smith #define __FUNCT__ "MGACycle_Private"
17dfbe8321SBarry Smith PetscErrorCode MGACycle_Private(MG *mg)
184b9ad928SBarry Smith {
196849ba73SBarry Smith   PetscErrorCode ierr;
20*79416396SBarry Smith   PetscInt       i,l = mg[0]->levels;
214b9ad928SBarry Smith   PetscScalar    zero = 0.0;
224b9ad928SBarry Smith 
234b9ad928SBarry Smith   PetscFunctionBegin;
244b9ad928SBarry Smith   /* compute RHS on each level */
254b9ad928SBarry Smith   for (i=l-1; i>0; i--) {
264b9ad928SBarry Smith     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
274b9ad928SBarry Smith   }
284b9ad928SBarry Smith   /* solve seperately on each level */
294b9ad928SBarry Smith   for (i=0; i<l; i++) {
304b9ad928SBarry Smith     ierr = VecSet(&zero,mg[i]->x);CHKERRQ(ierr);
314b9ad928SBarry Smith     if (mg[i]->eventsolve) {ierr = PetscLogEventBegin(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
3223ce1328SBarry Smith     ierr = KSPSolve(mg[i]->smoothd,mg[i]->b,mg[i]->x);CHKERRQ(ierr);
334b9ad928SBarry Smith     if (mg[i]->eventsolve) {ierr = PetscLogEventEnd(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
344b9ad928SBarry Smith   }
354b9ad928SBarry Smith   for (i=1; i<l; i++) {
364b9ad928SBarry Smith     ierr = MatInterpolateAdd(mg[i]->interpolate,mg[i-1]->x,mg[i]->x,mg[i]->x);CHKERRQ(ierr);
374b9ad928SBarry Smith   }
384b9ad928SBarry Smith   PetscFunctionReturn(0);
394b9ad928SBarry Smith }
40