1 /* 2 Full multigrid using either additive or multiplicative V or W cycle 3 */ 4 #include "src/ksp/pc/impls/mg/mgimpl.h" 5 6 EXTERN PetscErrorCode PCMGMCycle_Private(PC_MG **,PetscTruth*); 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "PCMGFCycle_Private" 10 PetscErrorCode PCMGFCycle_Private(PC_MG **mg) 11 { 12 PetscErrorCode ierr; 13 PetscInt i,l = mg[0]->levels; 14 15 PetscFunctionBegin; 16 /* restrict the RHS through all levels to coarsest. */ 17 for (i=l-1; i>0; i--){ 18 if (mg[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 19 ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr); 20 if (mg[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 21 } 22 23 /* work our way up through the levels */ 24 ierr = VecSet(mg[0]->x,0.0);CHKERRQ(ierr); 25 for (i=0; i<l-1; i++) { 26 ierr = PCMGMCycle_Private(&mg[i],PETSC_NULL);CHKERRQ(ierr); 27 if (mg[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 28 ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr); 29 if (mg[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 30 } 31 ierr = PCMGMCycle_Private(&mg[l-1],PETSC_NULL);CHKERRQ(ierr); 32 PetscFunctionReturn(0); 33 } 34 35 #undef __FUNCT__ 36 #define __FUNCT__ "PCMGKCycle_Private" 37 PetscErrorCode PCMGKCycle_Private(PC_MG **mg) 38 { 39 PetscErrorCode ierr; 40 PetscInt i,l = mg[0]->levels; 41 42 PetscFunctionBegin; 43 /* restrict the RHS through all levels to coarsest. */ 44 for (i=l-1; i>0; i--){ 45 if (mg[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 46 ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr); 47 if (mg[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 48 } 49 50 /* work our way up through the levels */ 51 ierr = VecSet(mg[0]->x,0.0);CHKERRQ(ierr); 52 for (i=0; i<l-1; i++) { 53 if (mg[i]->eventsmoothsolve) {ierr = PetscLogEventBegin(mg[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 54 ierr = KSPSolve(mg[i]->smoothd,mg[i]->b,mg[i]->x);CHKERRQ(ierr); 55 if (mg[i]->eventsmoothsolve) {ierr = PetscLogEventEnd(mg[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 56 if (mg[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 57 ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr); 58 if (mg[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mg[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 59 } 60 if (mg[l-1]->eventsmoothsolve) {ierr = PetscLogEventBegin(mg[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 61 ierr = KSPSolve(mg[l-1]->smoothd,mg[l-1]->b,mg[l-1]->x);CHKERRQ(ierr); 62 if (mg[l-1]->eventsmoothsolve) {ierr = PetscLogEventEnd(mg[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 63 64 PetscFunctionReturn(0); 65 } 66 67 68