xref: /petsc/src/ksp/pc/impls/mg/fmg.c (revision 31567311ea1f59b4e8129a5fbde7b11d292b5718)
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 {
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   /* restrict the RHS through all levels to coarsest. */
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   }
234b9ad928SBarry Smith 
244b9ad928SBarry Smith   /* work our way up through the levels */
25f3fbd535SBarry Smith   ierr = VecSet(mglevels[0]->x,0.0);CHKERRQ(ierr);
264b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
27f3fbd535SBarry Smith     ierr = PCMGMCycle_Private(pc,&mglevels[i],PETSC_NULL);CHKERRQ(ierr);
28*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
29f3fbd535SBarry Smith     ierr = MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);CHKERRQ(ierr);
30*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
314b9ad928SBarry Smith   }
32f3fbd535SBarry Smith   ierr = PCMGMCycle_Private(pc,&mglevels[l-1],PETSC_NULL);CHKERRQ(ierr);
334b9ad928SBarry Smith   PetscFunctionReturn(0);
344b9ad928SBarry Smith }
354b9ad928SBarry Smith 
364b9ad928SBarry Smith #undef __FUNCT__
379dcbbd2bSBarry Smith #define __FUNCT__ "PCMGKCycle_Private"
38*31567311SBarry Smith PetscErrorCode PCMGKCycle_Private(PC pc,PC_MG_Levels **mglevels)
394b9ad928SBarry Smith {
40*31567311SBarry Smith   PC_MG          *mg = (PC_MG*)pc->data;
416849ba73SBarry Smith   PetscErrorCode ierr;
42f3fbd535SBarry Smith   PetscInt       i,l = mglevels[0]->levels;
434b9ad928SBarry Smith 
444b9ad928SBarry Smith   PetscFunctionBegin;
454b9ad928SBarry Smith   /* restrict the RHS through all levels to coarsest. */
464b9ad928SBarry Smith   for (i=l-1; i>0; i--){
47*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
48f3fbd535SBarry Smith     ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr);
49*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
504b9ad928SBarry Smith   }
514b9ad928SBarry Smith 
524b9ad928SBarry Smith   /* work our way up through the levels */
53f3fbd535SBarry Smith   ierr = VecSet(mglevels[0]->x,0.0);CHKERRQ(ierr);
544b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
55*31567311SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
56f3fbd535SBarry Smith     ierr = KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);CHKERRQ(ierr);
57*31567311SBarry Smith     if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
58*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
59f3fbd535SBarry Smith     ierr = MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);CHKERRQ(ierr);
60*31567311SBarry Smith     if (mg->eventinterprestrict) {ierr = PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
614b9ad928SBarry Smith   }
62*31567311SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
63f3fbd535SBarry Smith   ierr = KSPSolve(mglevels[l-1]->smoothd,mglevels[l-1]->b,mglevels[l-1]->x);CHKERRQ(ierr);
64*31567311SBarry Smith   if (mg->eventsmoothsolve) {ierr = PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
654b9ad928SBarry Smith 
664b9ad928SBarry Smith   PetscFunctionReturn(0);
674b9ad928SBarry Smith }
684b9ad928SBarry Smith 
694b9ad928SBarry Smith 
70