xref: /petsc/src/ksp/pc/impls/mg/smg.c (revision 23ce13286abb4b24d31095dd761f3599be08cebe)
14b9ad928SBarry Smith /*$Id: smg.c,v 1.24 2001/08/07 03:03:36 balay Exp $*/
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith      Additive Multigrid V Cycle routine
44b9ad928SBarry Smith */
54b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h"
64b9ad928SBarry Smith 
74b9ad928SBarry Smith /*
84b9ad928SBarry Smith        MGACycle_Private - Given an MG structure created with MGCreate() runs
94b9ad928SBarry Smith                   one cycle down through the levels and back up. Applys
104b9ad928SBarry Smith                   the smoothers in an additive manner.
114b9ad928SBarry Smith 
124b9ad928SBarry Smith     Iput Parameters:
134b9ad928SBarry Smith .   mg - structure created with  MGCreate().
144b9ad928SBarry Smith 
154b9ad928SBarry Smith */
164b9ad928SBarry Smith #undef __FUNCT__
174b9ad928SBarry Smith #define __FUNCT__ "MGACycle_Private"
184b9ad928SBarry Smith int MGACycle_Private(MG *mg)
194b9ad928SBarry Smith {
204b9ad928SBarry Smith   int    i,l = mg[0]->levels,ierr;
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);}
32*23ce1328SBarry 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