xref: /petsc/src/ksp/pc/impls/mg/fmg.c (revision efb308899f58eff35e31f4136d1431f0c81bf4bc)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
44b9ad928SBarry Smith      Full multigrid using either additive or multiplicative V or W cycle
54b9ad928SBarry Smith */
64b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h"
74b9ad928SBarry Smith 
89dcbbd2bSBarry Smith EXTERN PetscErrorCode PCMGMCycle_Private(PC_MG **,PetscTruth*);
94b9ad928SBarry Smith 
104b9ad928SBarry Smith #undef __FUNCT__
119dcbbd2bSBarry Smith #define __FUNCT__ "PCMGFCycle_Private"
129dcbbd2bSBarry Smith PetscErrorCode PCMGFCycle_Private(PC_MG **mg)
134b9ad928SBarry Smith {
146849ba73SBarry Smith   PetscErrorCode ierr;
1579416396SBarry Smith   PetscInt       i,l = mg[0]->levels;
164b9ad928SBarry Smith 
174b9ad928SBarry Smith   PetscFunctionBegin;
184b9ad928SBarry Smith   /* restrict the RHS through all levels to coarsest. */
194b9ad928SBarry Smith   for (i=l-1; i>0; i--){
204b9ad928SBarry Smith     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
214b9ad928SBarry Smith   }
224b9ad928SBarry Smith 
234b9ad928SBarry Smith   /* work our way up through the levels */
24*efb30889SBarry Smith   ierr = VecSet(mg[0]->x,0.0);CHKERRQ(ierr);
254b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
269dcbbd2bSBarry Smith     ierr = PCMGMCycle_Private(&mg[i],PETSC_NULL);CHKERRQ(ierr);
274b9ad928SBarry Smith     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
284b9ad928SBarry Smith   }
299dcbbd2bSBarry Smith   ierr = PCMGMCycle_Private(&mg[l-1],PETSC_NULL);CHKERRQ(ierr);
304b9ad928SBarry Smith   PetscFunctionReturn(0);
314b9ad928SBarry Smith }
324b9ad928SBarry Smith 
334b9ad928SBarry Smith #undef __FUNCT__
349dcbbd2bSBarry Smith #define __FUNCT__ "PCMGKCycle_Private"
359dcbbd2bSBarry Smith PetscErrorCode PCMGKCycle_Private(PC_MG **mg)
364b9ad928SBarry Smith {
376849ba73SBarry Smith   PetscErrorCode ierr;
3879416396SBarry Smith   PetscInt       i,l = mg[0]->levels;
394b9ad928SBarry Smith 
404b9ad928SBarry Smith   PetscFunctionBegin;
414b9ad928SBarry Smith   /* restrict the RHS through all levels to coarsest. */
424b9ad928SBarry Smith   for (i=l-1; i>0; i--){
434b9ad928SBarry Smith     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
444b9ad928SBarry Smith   }
454b9ad928SBarry Smith 
464b9ad928SBarry Smith   /* work our way up through the levels */
47*efb30889SBarry Smith   ierr = VecSet(mg[0]->x,0.0);CHKERRQ(ierr);
484b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
494b9ad928SBarry Smith     if (mg[i]->eventsolve) {ierr = PetscLogEventBegin(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
5023ce1328SBarry Smith     ierr = KSPSolve(mg[i]->smoothd,mg[i]->b,mg[i]->x);CHKERRQ(ierr);
514b9ad928SBarry Smith     if (mg[i]->eventsolve) {ierr = PetscLogEventEnd(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
524b9ad928SBarry Smith     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
534b9ad928SBarry Smith   }
544b9ad928SBarry Smith   if (mg[l-1]->eventsolve) {ierr = PetscLogEventBegin(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
5523ce1328SBarry Smith   ierr = KSPSolve(mg[l-1]->smoothd,mg[l-1]->b,mg[l-1]->x);CHKERRQ(ierr);
564b9ad928SBarry Smith   if (mg[l-1]->eventsolve) {ierr = PetscLogEventEnd(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith 
584b9ad928SBarry Smith   PetscFunctionReturn(0);
594b9ad928SBarry Smith }
604b9ad928SBarry Smith 
614b9ad928SBarry Smith 
62