Actual source code: fmg.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: /*
  2:      Full multigrid using either additive or multiplicative V or W cycle
  3: */
  4:  #include <petsc/private/pcmgimpl.h>

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

  8: PetscErrorCode PCMGFCycle_Private(PC pc,PC_MG_Levels **mglevels)
  9: {
 11:   PetscInt       i,l = mglevels[0]->levels;

 14:   /* restrict the RHS through all levels to coarsest. */
 15:   for (i=l-1; i>0; i--) {
 16:     if (mglevels[i]->eventinterprestrict) {PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);}
 17:     MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);
 18:     if (mglevels[i]->eventinterprestrict) {PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);}
 19:   }

 21:   /* work our way up through the levels */
 22:   VecSet(mglevels[0]->x,0.0);
 23:   for (i=0; i<l-1; i++) {
 24:     PCMGMCycle_Private(pc,&mglevels[i],NULL);
 25:     if (mglevels[i+1]->eventinterprestrict) {PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);}
 26:     MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);
 27:     if (mglevels[i+1]->eventinterprestrict) {PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);}
 28:   }
 29:   PCMGMCycle_Private(pc,&mglevels[l-1],NULL);
 30:   return(0);
 31: }

 33: PetscErrorCode PCMGKCycle_Private(PC pc,PC_MG_Levels **mglevels)
 34: {
 36:   PetscInt       i,l = mglevels[0]->levels;

 39:   /* restrict the RHS through all levels to coarsest. */
 40:   for (i=l-1; i>0; i--) {
 41:     if (mglevels[i]->eventinterprestrict) {PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);}
 42:     MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);
 43:     if (mglevels[i]->eventinterprestrict) {PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);}
 44:   }

 46:   /* work our way up through the levels */
 47:   VecSet(mglevels[0]->x,0.0);
 48:   for (i=0; i<l-1; i++) {
 49:     if (mglevels[i]->eventsmoothsolve) {PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0);}
 50:     KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);
 51:     KSPCheckSolve(mglevels[i]->smoothd,pc,mglevels[i]->x);
 52:     if (mglevels[i]->eventsmoothsolve) {PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0);}
 53:     if (mglevels[i+1]->eventinterprestrict) {PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);}
 54:     MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);
 55:     if (mglevels[i+1]->eventinterprestrict) {PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);}
 56:   }
 57:   if (mglevels[l-1]->eventsmoothsolve) {PetscLogEventBegin(mglevels[l-1]->eventsmoothsolve,0,0,0,0);}
 58:   KSPSolve(mglevels[l-1]->smoothd,mglevels[l-1]->b,mglevels[l-1]->x);
 59:   KSPCheckSolve(mglevels[l-1]->smoothd,pc,mglevels[l-1]->x);
 60:   if (mglevels[l-1]->eventsmoothsolve) {PetscLogEventEnd(mglevels[l-1]->eventsmoothsolve,0,0,0,0);}
 61:   return(0);
 62: }