xref: /petsc/src/ksp/pc/impls/mg/fmg.c (revision fcb023d4cd21af6e7bfc4dfb4b9546d75029f707)
14b9ad928SBarry Smith /*
24b9ad928SBarry Smith      Full multigrid using either additive or multiplicative V or W cycle
34b9ad928SBarry Smith */
4af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>
54b9ad928SBarry Smith 
6*fcb023d4SJed Brown PetscErrorCode PCMGFCycle_Private(PC pc,PC_MG_Levels **mglevels,PetscBool transpose)
74b9ad928SBarry Smith {
86849ba73SBarry Smith   PetscErrorCode ierr;
9f3fbd535SBarry Smith   PetscInt       i,l = mglevels[0]->levels;
104b9ad928SBarry Smith 
114b9ad928SBarry Smith   PetscFunctionBegin;
12*fcb023d4SJed Brown   if (!transpose) {
134b9ad928SBarry Smith     /* restrict the RHS through all levels to coarsest. */
144b9ad928SBarry Smith     for (i=l-1; i>0; i--) {
1563e6d426SJed Brown       if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
16f3fbd535SBarry Smith       ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr);
1763e6d426SJed Brown       if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
184b9ad928SBarry Smith     }
194b9ad928SBarry Smith 
204b9ad928SBarry Smith     /* work our way up through the levels */
21f3fbd535SBarry Smith     ierr = VecSet(mglevels[0]->x,0.0);CHKERRQ(ierr);
224b9ad928SBarry Smith     for (i=0; i<l-1; i++) {
23*fcb023d4SJed Brown       ierr = PCMGMCycle_Private(pc,&mglevels[i],transpose,NULL);CHKERRQ(ierr);
2463e6d426SJed Brown       if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
25f3fbd535SBarry Smith       ierr = MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);CHKERRQ(ierr);
2663e6d426SJed Brown       if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
274b9ad928SBarry Smith     }
28*fcb023d4SJed Brown     ierr = PCMGMCycle_Private(pc,&mglevels[l-1],transpose,NULL);CHKERRQ(ierr);
29*fcb023d4SJed Brown   } else {
30*fcb023d4SJed Brown     ierr = PCMGMCycle_Private(pc,&mglevels[l-1],transpose,NULL);CHKERRQ(ierr);
31*fcb023d4SJed Brown     for (i=l-2; i>=0; i--) {
32*fcb023d4SJed Brown       if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
33*fcb023d4SJed Brown       ierr = MatRestrict(mglevels[i+1]->interpolate,mglevels[i+1]->x,mglevels[i]->x);CHKERRQ(ierr);
34*fcb023d4SJed Brown       if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
35*fcb023d4SJed Brown       ierr = PCMGMCycle_Private(pc,&mglevels[i],transpose,NULL);CHKERRQ(ierr);
36*fcb023d4SJed Brown     }
37*fcb023d4SJed Brown     for (i=1; i<l; i++) {
38*fcb023d4SJed Brown       if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
39*fcb023d4SJed Brown       ierr = MatInterpolate(mglevels[i]->restrct,mglevels[i-1]->b,mglevels[i]->b);CHKERRQ(ierr);
40*fcb023d4SJed Brown       if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
41*fcb023d4SJed Brown     }
42*fcb023d4SJed Brown   }
434b9ad928SBarry Smith   PetscFunctionReturn(0);
444b9ad928SBarry Smith }
454b9ad928SBarry Smith 
46*fcb023d4SJed Brown PetscErrorCode PCMGKCycle_Private(PC pc,PC_MG_Levels **mglevels,PetscBool transpose)
474b9ad928SBarry Smith {
486849ba73SBarry Smith   PetscErrorCode ierr;
49f3fbd535SBarry Smith   PetscInt       i,l = mglevels[0]->levels;
504b9ad928SBarry Smith 
514b9ad928SBarry Smith   PetscFunctionBegin;
524b9ad928SBarry Smith   /* restrict the RHS through all levels to coarsest. */
534b9ad928SBarry Smith   for (i=l-1; i>0; i--) {
5463e6d426SJed Brown     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
55f3fbd535SBarry Smith     ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr);
5663e6d426SJed Brown     if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
574b9ad928SBarry Smith   }
584b9ad928SBarry Smith 
594b9ad928SBarry Smith   /* work our way up through the levels */
60f3fbd535SBarry Smith   ierr = VecSet(mglevels[0]->x,0.0);CHKERRQ(ierr);
614b9ad928SBarry Smith   for (i=0; i<l-1; i++) {
6263e6d426SJed Brown     if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
63f3fbd535SBarry Smith     ierr = KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);CHKERRQ(ierr);
64c0decd05SBarry Smith     ierr = KSPCheckSolve(mglevels[i]->smoothd,pc,mglevels[i]->x);CHKERRQ(ierr);
6563e6d426SJed Brown     if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
6663e6d426SJed Brown     if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
67f3fbd535SBarry Smith     ierr = MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);CHKERRQ(ierr);
6863e6d426SJed Brown     if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
694b9ad928SBarry Smith   }
7063e6d426SJed Brown   if (mglevels[l-1]->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
71f3fbd535SBarry Smith   ierr = KSPSolve(mglevels[l-1]->smoothd,mglevels[l-1]->b,mglevels[l-1]->x);CHKERRQ(ierr);
72c0decd05SBarry Smith   ierr = KSPCheckSolve(mglevels[l-1]->smoothd,pc,mglevels[l-1]->x);CHKERRQ(ierr);
7363e6d426SJed Brown   if (mglevels[l-1]->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
744b9ad928SBarry Smith   PetscFunctionReturn(0);
754b9ad928SBarry Smith }
764b9ad928SBarry Smith 
774b9ad928SBarry Smith 
78