1 #define PETSCKSP_DLL 2 3 /* 4 Additive Multigrid V Cycle routine 5 */ 6 #include "src/ksp/pc/impls/mg/mgimpl.h" 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "PCMGACycle_Private" 10 PetscErrorCode PCMGACycle_Private(PC_MG **mg) 11 { 12 PetscErrorCode ierr; 13 PetscInt i,l = mg[0]->levels; 14 15 PetscFunctionBegin; 16 /* compute RHS on each level */ 17 for (i=l-1; i>0; i--) { 18 ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr); 19 } 20 /* solve seperately on each level */ 21 for (i=0; i<l; i++) { 22 ierr = VecSet(mg[i]->x,0.0);CHKERRQ(ierr); 23 if (mg[i]->eventsolve) {ierr = PetscLogEventBegin(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);} 24 ierr = KSPSolve(mg[i]->smoothd,mg[i]->b,mg[i]->x);CHKERRQ(ierr); 25 if (mg[i]->eventsolve) {ierr = PetscLogEventEnd(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);} 26 } 27 for (i=1; i<l; i++) { 28 ierr = MatInterpolateAdd(mg[i]->interpolate,mg[i-1]->x,mg[i]->x,mg[i]->x);CHKERRQ(ierr); 29 } 30 PetscFunctionReturn(0); 31 } 32