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