xref: /petsc/src/ksp/pc/impls/mg/smg.c (revision 31567311ea1f59b4e8129a5fbde7b11d292b5718)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith      Additive Multigrid V Cycle routine
54b9ad928SBarry Smith */
67c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #undef __FUNCT__
99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGACycle_Private"
10*31567311SBarry Smith PetscErrorCode PCMGACycle_Private(PC pc,PC_MG_Levels **mglevels)
114b9ad928SBarry Smith {
12*31567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
136849ba73SBarry Smith   PetscErrorCode ierr;
14f3fbd535SBarry Smith   PetscInt       i,l = mglevels[0]->levels;
154b9ad928SBarry Smith 
164b9ad928SBarry Smith   PetscFunctionBegin;
174b9ad928SBarry Smith   /* compute RHS on each level */
184b9ad928SBarry Smith   for (i=l-1; i>0; i--) {
19*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
20f3fbd535SBarry Smith     ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr);
21*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
224b9ad928SBarry Smith   }
23a8c7a070SBarry Smith   /* solve separately on each level */
244b9ad928SBarry Smith   for (i=0; i<l; i++) {
25f3fbd535SBarry Smith     ierr = VecSet(mglevels[i]->x,0.0);CHKERRQ(ierr);
26*31567311SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
27f3fbd535SBarry Smith     ierr = KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);CHKERRQ(ierr);
28*31567311SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
294b9ad928SBarry Smith   }
304b9ad928SBarry Smith   for (i=1; i<l; i++) {
31*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
32f3fbd535SBarry Smith     ierr = MatInterpolateAdd(mglevels[i]->interpolate,mglevels[i-1]->x,mglevels[i]->x,mglevels[i]->x);CHKERRQ(ierr);
33*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
344b9ad928SBarry Smith   }
354b9ad928SBarry Smith   PetscFunctionReturn(0);
364b9ad928SBarry Smith }
37