xref: /petsc/src/ksp/pc/impls/mg/fmg.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith      Full multigrid using either additive or multiplicative V or W cycle
34b9ad928SBarry Smith */
44b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h"
54b9ad928SBarry Smith 
6dfbe8321SBarry Smith EXTERN PetscErrorCode MGMCycle_Private(MG *,PetscTruth*);
74b9ad928SBarry Smith 
84b9ad928SBarry Smith /*
94b9ad928SBarry Smith        MGFCycle_Private - Given an MG structure created with MGCreate() runs
104b9ad928SBarry Smith                full multigrid.
114b9ad928SBarry Smith 
124b9ad928SBarry Smith     Iput Parameters:
134b9ad928SBarry Smith .   mg - structure created with MGCreate().
144b9ad928SBarry Smith 
154b9ad928SBarry Smith     Note: This may not be what others call full multigrid. What we
164b9ad928SBarry Smith           do is restrict the rhs to all levels, then starting
174b9ad928SBarry Smith           on the coarsest level work our way up generating
184b9ad928SBarry Smith           initial guess for the next level. This provides an
194b9ad928SBarry Smith           improved preconditioner but not a great improvement.
204b9ad928SBarry Smith */
214b9ad928SBarry Smith #undef __FUNCT__
224b9ad928SBarry Smith #define __FUNCT__ "MGFCycle_Private"
23dfbe8321SBarry Smith PetscErrorCode MGFCycle_Private(MG *mg)
244b9ad928SBarry Smith {
25*6849ba73SBarry Smith   PetscErrorCode ierr;
26*6849ba73SBarry Smith   int    i,l = mg[0]->levels;
274b9ad928SBarry Smith   PetscScalar zero = 0.0;
284b9ad928SBarry Smith 
294b9ad928SBarry Smith   PetscFunctionBegin;
304b9ad928SBarry Smith   /* restrict the RHS through all levels to coarsest. */
314b9ad928SBarry Smith   for (i=l-1; i>0; i--){
324b9ad928SBarry Smith     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
334b9ad928SBarry Smith   }
344b9ad928SBarry Smith 
354b9ad928SBarry Smith   /* work our way up through the levels */
364b9ad928SBarry Smith   ierr = VecSet(&zero,mg[0]->x);CHKERRQ(ierr);
374b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
384b9ad928SBarry Smith     ierr = MGMCycle_Private(&mg[i],PETSC_NULL);CHKERRQ(ierr);
394b9ad928SBarry Smith     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
404b9ad928SBarry Smith   }
414b9ad928SBarry Smith   ierr = MGMCycle_Private(&mg[l-1],PETSC_NULL);CHKERRQ(ierr);
424b9ad928SBarry Smith   PetscFunctionReturn(0);
434b9ad928SBarry Smith }
444b9ad928SBarry Smith 
454b9ad928SBarry Smith /*
464b9ad928SBarry Smith        MGKCycle_Private - Given an MG structure created with MGCreate() runs
474b9ad928SBarry Smith                full Kascade MG solve.
484b9ad928SBarry Smith 
494b9ad928SBarry Smith     Iput Parameters:
504b9ad928SBarry Smith .   mg - structure created with MGCreate().
514b9ad928SBarry Smith 
524b9ad928SBarry Smith     Note: This may not be what others call Kascadic MG.
534b9ad928SBarry Smith */
544b9ad928SBarry Smith #undef __FUNCT__
554b9ad928SBarry Smith #define __FUNCT__ "MGKCycle_Private"
56dfbe8321SBarry Smith PetscErrorCode MGKCycle_Private(MG *mg)
574b9ad928SBarry Smith {
58*6849ba73SBarry Smith   PetscErrorCode ierr;
59*6849ba73SBarry Smith   int    i,l = mg[0]->levels;
604b9ad928SBarry Smith   PetscScalar zero = 0.0;
614b9ad928SBarry Smith 
624b9ad928SBarry Smith   PetscFunctionBegin;
634b9ad928SBarry Smith   /* restrict the RHS through all levels to coarsest. */
644b9ad928SBarry Smith   for (i=l-1; i>0; i--){
654b9ad928SBarry Smith     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
664b9ad928SBarry Smith   }
674b9ad928SBarry Smith 
684b9ad928SBarry Smith   /* work our way up through the levels */
694b9ad928SBarry Smith   ierr = VecSet(&zero,mg[0]->x);CHKERRQ(ierr);
704b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
714b9ad928SBarry Smith     if (mg[i]->eventsolve) {ierr = PetscLogEventBegin(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
7223ce1328SBarry Smith     ierr = KSPSolve(mg[i]->smoothd,mg[i]->b,mg[i]->x);CHKERRQ(ierr);
734b9ad928SBarry Smith     if (mg[i]->eventsolve) {ierr = PetscLogEventEnd(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
744b9ad928SBarry Smith     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
754b9ad928SBarry Smith   }
764b9ad928SBarry Smith   if (mg[l-1]->eventsolve) {ierr = PetscLogEventBegin(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
7723ce1328SBarry Smith   ierr = KSPSolve(mg[l-1]->smoothd,mg[l-1]->b,mg[l-1]->x);CHKERRQ(ierr);
784b9ad928SBarry Smith   if (mg[l-1]->eventsolve) {ierr = PetscLogEventEnd(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
794b9ad928SBarry Smith 
804b9ad928SBarry Smith   PetscFunctionReturn(0);
814b9ad928SBarry Smith }
824b9ad928SBarry Smith 
834b9ad928SBarry Smith 
84