xref: /petsc/src/ksp/pc/impls/mg/fmg.c (revision 7c4f633dc6bb6149cca88d301ead35a99e103cbb)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith      Full multigrid using either additive or multiplicative V or W cycle
34b9ad928SBarry Smith */
4*7c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"
54b9ad928SBarry Smith 
6b7d41711SBarry Smith EXTERN PetscErrorCode PCMGMCycle_Private(PC,PC_MG **,PCRichardsonConvergedReason*);
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #undef __FUNCT__
99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGFCycle_Private"
101e2582c4SBarry Smith PetscErrorCode PCMGFCycle_Private(PC pc,PC_MG **mg)
114b9ad928SBarry Smith {
126849ba73SBarry Smith   PetscErrorCode ierr;
1379416396SBarry Smith   PetscInt       i,l = mg[0]->levels;
144b9ad928SBarry Smith 
154b9ad928SBarry Smith   PetscFunctionBegin;
164b9ad928SBarry Smith   /* restrict the RHS through all levels to coarsest. */
174b9ad928SBarry Smith   for (i=l-1; i>0; i--){
1832cf1786SBarry Smith     if (mg[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
194b9ad928SBarry Smith     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
2032cf1786SBarry Smith     if (mg[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
214b9ad928SBarry Smith   }
224b9ad928SBarry Smith 
234b9ad928SBarry Smith   /* work our way up through the levels */
24efb30889SBarry Smith   ierr = VecSet(mg[0]->x,0.0);CHKERRQ(ierr);
254b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
261e2582c4SBarry Smith     ierr = PCMGMCycle_Private(pc,&mg[i],PETSC_NULL);CHKERRQ(ierr);
2732cf1786SBarry Smith     if (mg[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
284b9ad928SBarry Smith     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
2932cf1786SBarry Smith     if (mg[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
304b9ad928SBarry Smith   }
311e2582c4SBarry Smith   ierr = PCMGMCycle_Private(pc,&mg[l-1],PETSC_NULL);CHKERRQ(ierr);
324b9ad928SBarry Smith   PetscFunctionReturn(0);
334b9ad928SBarry Smith }
344b9ad928SBarry Smith 
354b9ad928SBarry Smith #undef __FUNCT__
369dcbbd2bSBarry Smith #define __FUNCT__ "PCMGKCycle_Private"
379dcbbd2bSBarry Smith PetscErrorCode PCMGKCycle_Private(PC_MG **mg)
384b9ad928SBarry Smith {
396849ba73SBarry Smith   PetscErrorCode ierr;
4079416396SBarry Smith   PetscInt       i,l = mg[0]->levels;
414b9ad928SBarry Smith 
424b9ad928SBarry Smith   PetscFunctionBegin;
434b9ad928SBarry Smith   /* restrict the RHS through all levels to coarsest. */
444b9ad928SBarry Smith   for (i=l-1; i>0; i--){
4532cf1786SBarry Smith     if (mg[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
464b9ad928SBarry Smith     ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr);
4732cf1786SBarry Smith     if (mg[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
484b9ad928SBarry Smith   }
494b9ad928SBarry Smith 
504b9ad928SBarry Smith   /* work our way up through the levels */
51efb30889SBarry Smith   ierr = VecSet(mg[0]->x,0.0);CHKERRQ(ierr);
524b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
5332cf1786SBarry Smith     if (mg[i]->eventsmoothsolve) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5423ce1328SBarry Smith     ierr = KSPSolve(mg[i]->smoothd,mg[i]->b,mg[i]->x);CHKERRQ(ierr);
5532cf1786SBarry Smith     if (mg[i]->eventsmoothsolve) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
5632cf1786SBarry Smith     if (mg[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith     ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr);
5832cf1786SBarry Smith     if (mg[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
594b9ad928SBarry Smith   }
6032cf1786SBarry Smith   if (mg[l-1]->eventsmoothsolve) {ierr = PetscLogEventBegin(mg[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
6123ce1328SBarry Smith   ierr = KSPSolve(mg[l-1]->smoothd,mg[l-1]->b,mg[l-1]->x);CHKERRQ(ierr);
6232cf1786SBarry Smith   if (mg[l-1]->eventsmoothsolve) {ierr = PetscLogEventEnd(mg[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
634b9ad928SBarry Smith 
644b9ad928SBarry Smith   PetscFunctionReturn(0);
654b9ad928SBarry Smith }
664b9ad928SBarry Smith 
674b9ad928SBarry Smith 
68