Actual source code: fmg.c

  1: /*
  2:      Full multigrid using either additive or multiplicative V or W cycle
  3: */
 4:  #include ../src/ksp/pc/impls/mg/mgimpl.h

  6: EXTERN PetscErrorCode PCMGMCycle_Private(PC,PC_MG_Levels **,PCRichardsonConvergedReason*);

 10: PetscErrorCode PCMGFCycle_Private(PC pc,PC_MG_Levels **mglevels)
 11: {
 12:   PC_MG          *mg = (PC_MG*)pc->data;
 14:   PetscInt       i,l = mglevels[0]->levels;

 17:   /* restrict the RHS through all levels to coarsest. */
 18:   for (i=l-1; i>0; i--){
 19:     if (mg->eventinterprestrict) {PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);}
 20:     MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);
 21:     if (mg->eventinterprestrict) {PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);}
 22:   }
 23: 
 24:   /* work our way up through the levels */
 25:   VecSet(mglevels[0]->x,0.0);
 26:   for (i=0; i<l-1; i++) {
 27:     PCMGMCycle_Private(pc,&mglevels[i],PETSC_NULL);
 28:     if (mg->eventinterprestrict) {PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);}
 29:     MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);
 30:     if (mg->eventinterprestrict) {PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);}
 31:   }
 32:   PCMGMCycle_Private(pc,&mglevels[l-1],PETSC_NULL);
 33:   return(0);
 34: }

 38: PetscErrorCode PCMGKCycle_Private(PC pc,PC_MG_Levels **mglevels)
 39: {
 40:   PC_MG          *mg = (PC_MG*)pc->data;
 42:   PetscInt       i,l = mglevels[0]->levels;

 45:   /* restrict the RHS through all levels to coarsest. */
 46:   for (i=l-1; i>0; i--){
 47:     if (mg->eventinterprestrict) {PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);}
 48:     MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);
 49:     if (mg->eventinterprestrict) {PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);}
 50:   }
 51: 
 52:   /* work our way up through the levels */
 53:   VecSet(mglevels[0]->x,0.0);
 54:   for (i=0; i<l-1; i++) {
 55:     if (mg->eventsmoothsolve) {PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);}
 56:     KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);
 57:     if (mg->eventsmoothsolve) {PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);}
 58:     if (mg->eventinterprestrict) {PetscLogEventBegin(mg->eventinterprestrict,0,0,0,0);}
 59:     MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);
 60:     if (mg->eventinterprestrict) {PetscLogEventEnd(mg->eventinterprestrict,0,0,0,0);}
 61:   }
 62:   if (mg->eventsmoothsolve) {PetscLogEventBegin(mg->eventsmoothsolve,0,0,0,0);}
 63:   KSPSolve(mglevels[l-1]->smoothd,mglevels[l-1]->b,mglevels[l-1]->x);
 64:   if (mg->eventsmoothsolve) {PetscLogEventEnd(mg->eventsmoothsolve,0,0,0,0);}

 66:   return(0);
 67: }