Actual source code: smg.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2: /*
  3:      Additive Multigrid V Cycle routine
  4: */
  5:  #include <petsc/private/pcmgimpl.h>

  7: PetscErrorCode PCMGACycle_Private(PC pc,PC_MG_Levels **mglevels)
  8: {
 10:   PetscInt       i,l = mglevels[0]->levels;

 13:   /* compute RHS on each level */
 14:   for (i=l-1; i>0; i--) {
 15:     if (mglevels[i]->eventinterprestrict) {PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);}
 16:     MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);
 17:     if (mglevels[i]->eventinterprestrict) {PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);}
 18:   }
 19:   /* solve separately on each level */
 20:   for (i=0; i<l; i++) {
 21:     VecSet(mglevels[i]->x,0.0);
 22:     if (mglevels[i]->eventsmoothsolve) {PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0);}
 23:     KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);
 24:     KSPCheckSolve(mglevels[i]->smoothd,pc,mglevels[i]->x);
 25:     if (mglevels[i]->eventsmoothsolve) {PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0);}
 26:   }
 27:   for (i=1; i<l; i++) {
 28:     if (mglevels[i]->eventinterprestrict) {PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);}
 29:     MatInterpolateAdd(mglevels[i]->interpolate,mglevels[i-1]->x,mglevels[i]->x,mglevels[i]->x);
 30:     if (mglevels[i]->eventinterprestrict) {PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);}
 31:   }
 32:   return(0);
 33: }