xref: /petsc/src/ksp/pc/impls/mg/fmg.c (revision 63e6d426ec67d9b55d8f033309a2f8ed2bb643f7)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith      Full multigrid using either additive or multiplicative V or W cycle
34b9ad928SBarry Smith */
47c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h"
54b9ad928SBarry Smith 
6f3fbd535SBarry Smith EXTERN PetscErrorCode PCMGMCycle_Private(PC,PC_MG_Levels **,PCRichardsonConvergedReason*);
74b9ad928SBarry Smith 
84b9ad928SBarry Smith #undef __FUNCT__
99dcbbd2bSBarry Smith #define __FUNCT__ "PCMGFCycle_Private"
10f3fbd535SBarry Smith PetscErrorCode PCMGFCycle_Private(PC pc,PC_MG_Levels **mglevels)
114b9ad928SBarry Smith {
126849ba73SBarry Smith   PetscErrorCode ierr;
13f3fbd535SBarry Smith   PetscInt       i,l = mglevels[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--){
18*63e6d426SJed Brown     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
19f3fbd535SBarry Smith     ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr);
20*63e6d426SJed Brown     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
214b9ad928SBarry Smith   }
224b9ad928SBarry Smith 
234b9ad928SBarry Smith   /* work our way up through the levels */
24f3fbd535SBarry Smith   ierr = VecSet(mglevels[0]->x,0.0);CHKERRQ(ierr);
254b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
26f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,&mglevels[i],PETSC_NULL);CHKERRQ(ierr);
27*63e6d426SJed Brown     if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
28f3fbd535SBarry Smith     ierr = MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);CHKERRQ(ierr);
29*63e6d426SJed Brown     if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
304b9ad928SBarry Smith   }
31f3fbd535SBarry Smith   ierr = PCMGMCycle_Private(pc,&mglevels[l-1],PETSC_NULL);CHKERRQ(ierr);
324b9ad928SBarry Smith   PetscFunctionReturn(0);
334b9ad928SBarry Smith }
344b9ad928SBarry Smith 
354b9ad928SBarry Smith #undef __FUNCT__
369dcbbd2bSBarry Smith #define __FUNCT__ "PCMGKCycle_Private"
3731567311SBarry Smith PetscErrorCode PCMGKCycle_Private(PC pc,PC_MG_Levels **mglevels)
384b9ad928SBarry Smith {
396849ba73SBarry Smith   PetscErrorCode ierr;
40f3fbd535SBarry Smith   PetscInt       i,l = mglevels[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--){
45*63e6d426SJed Brown     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
46f3fbd535SBarry Smith     ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr);
47*63e6d426SJed Brown     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
484b9ad928SBarry Smith   }
494b9ad928SBarry Smith 
504b9ad928SBarry Smith   /* work our way up through the levels */
51f3fbd535SBarry Smith   ierr = VecSet(mglevels[0]->x,0.0);CHKERRQ(ierr);
524b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
53*63e6d426SJed Brown     if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
54f3fbd535SBarry Smith     ierr = KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);CHKERRQ(ierr);
55*63e6d426SJed Brown     if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
56*63e6d426SJed Brown     if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
57f3fbd535SBarry Smith     ierr = MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);CHKERRQ(ierr);
58*63e6d426SJed Brown     if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
594b9ad928SBarry Smith   }
60*63e6d426SJed Brown   if (mglevels[l-1]->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
61f3fbd535SBarry Smith   ierr = KSPSolve(mglevels[l-1]->smoothd,mglevels[l-1]->b,mglevels[l-1]->x);CHKERRQ(ierr);
62*63e6d426SJed Brown   if (mglevels[l-1]->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
634b9ad928SBarry Smith 
644b9ad928SBarry Smith   PetscFunctionReturn(0);
654b9ad928SBarry Smith }
664b9ad928SBarry Smith 
674b9ad928SBarry Smith 
68