1*4b9ad928SBarry Smith /*$Id: fmg.c,v 1.26 2001/08/21 21:03:20 bsmith Exp $*/ 2*4b9ad928SBarry Smith /* 3*4b9ad928SBarry Smith Full multigrid using either additive or multiplicative V or W cycle 4*4b9ad928SBarry Smith */ 5*4b9ad928SBarry Smith #include "src/ksp/pc/impls/mg/mgimpl.h" 6*4b9ad928SBarry Smith 7*4b9ad928SBarry Smith EXTERN int MGMCycle_Private(MG *,PetscTruth*); 8*4b9ad928SBarry Smith 9*4b9ad928SBarry Smith /* 10*4b9ad928SBarry Smith MGFCycle_Private - Given an MG structure created with MGCreate() runs 11*4b9ad928SBarry Smith full multigrid. 12*4b9ad928SBarry Smith 13*4b9ad928SBarry Smith Iput Parameters: 14*4b9ad928SBarry Smith . mg - structure created with MGCreate(). 15*4b9ad928SBarry Smith 16*4b9ad928SBarry Smith Note: This may not be what others call full multigrid. What we 17*4b9ad928SBarry Smith do is restrict the rhs to all levels, then starting 18*4b9ad928SBarry Smith on the coarsest level work our way up generating 19*4b9ad928SBarry Smith initial guess for the next level. This provides an 20*4b9ad928SBarry Smith improved preconditioner but not a great improvement. 21*4b9ad928SBarry Smith */ 22*4b9ad928SBarry Smith #undef __FUNCT__ 23*4b9ad928SBarry Smith #define __FUNCT__ "MGFCycle_Private" 24*4b9ad928SBarry Smith int MGFCycle_Private(MG *mg) 25*4b9ad928SBarry Smith { 26*4b9ad928SBarry Smith int i,l = mg[0]->levels,ierr; 27*4b9ad928SBarry Smith PetscScalar zero = 0.0; 28*4b9ad928SBarry Smith 29*4b9ad928SBarry Smith PetscFunctionBegin; 30*4b9ad928SBarry Smith /* restrict the RHS through all levels to coarsest. */ 31*4b9ad928SBarry Smith for (i=l-1; i>0; i--){ 32*4b9ad928SBarry Smith ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr); 33*4b9ad928SBarry Smith } 34*4b9ad928SBarry Smith 35*4b9ad928SBarry Smith /* work our way up through the levels */ 36*4b9ad928SBarry Smith ierr = VecSet(&zero,mg[0]->x);CHKERRQ(ierr); 37*4b9ad928SBarry Smith for (i=0; i<l-1; i++) { 38*4b9ad928SBarry Smith ierr = MGMCycle_Private(&mg[i],PETSC_NULL);CHKERRQ(ierr); 39*4b9ad928SBarry Smith ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr); 40*4b9ad928SBarry Smith } 41*4b9ad928SBarry Smith ierr = MGMCycle_Private(&mg[l-1],PETSC_NULL);CHKERRQ(ierr); 42*4b9ad928SBarry Smith PetscFunctionReturn(0); 43*4b9ad928SBarry Smith } 44*4b9ad928SBarry Smith 45*4b9ad928SBarry Smith /* 46*4b9ad928SBarry Smith MGKCycle_Private - Given an MG structure created with MGCreate() runs 47*4b9ad928SBarry Smith full Kascade MG solve. 48*4b9ad928SBarry Smith 49*4b9ad928SBarry Smith Iput Parameters: 50*4b9ad928SBarry Smith . mg - structure created with MGCreate(). 51*4b9ad928SBarry Smith 52*4b9ad928SBarry Smith Note: This may not be what others call Kascadic MG. 53*4b9ad928SBarry Smith */ 54*4b9ad928SBarry Smith #undef __FUNCT__ 55*4b9ad928SBarry Smith #define __FUNCT__ "MGKCycle_Private" 56*4b9ad928SBarry Smith int MGKCycle_Private(MG *mg) 57*4b9ad928SBarry Smith { 58*4b9ad928SBarry Smith int i,l = mg[0]->levels,ierr; 59*4b9ad928SBarry Smith PetscScalar zero = 0.0; 60*4b9ad928SBarry Smith 61*4b9ad928SBarry Smith PetscFunctionBegin; 62*4b9ad928SBarry Smith /* restrict the RHS through all levels to coarsest. */ 63*4b9ad928SBarry Smith for (i=l-1; i>0; i--){ 64*4b9ad928SBarry Smith ierr = MatRestrict(mg[i]->restrct,mg[i]->b,mg[i-1]->b);CHKERRQ(ierr); 65*4b9ad928SBarry Smith } 66*4b9ad928SBarry Smith 67*4b9ad928SBarry Smith /* work our way up through the levels */ 68*4b9ad928SBarry Smith ierr = VecSet(&zero,mg[0]->x);CHKERRQ(ierr); 69*4b9ad928SBarry Smith for (i=0; i<l-1; i++) { 70*4b9ad928SBarry Smith if (mg[i]->eventsolve) {ierr = PetscLogEventBegin(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);} 71*4b9ad928SBarry Smith ierr = KSPSetRhs(mg[i]->smoothd,mg[i]->b);CHKERRQ(ierr); 72*4b9ad928SBarry Smith ierr = KSPSetSolution(mg[i]->smoothd,mg[i]->x);CHKERRQ(ierr); 73*4b9ad928SBarry Smith ierr = KSPSolve(mg[i]->smoothd);CHKERRQ(ierr); 74*4b9ad928SBarry Smith if (mg[i]->eventsolve) {ierr = PetscLogEventEnd(mg[i]->eventsolve,0,0,0,0);CHKERRQ(ierr);} 75*4b9ad928SBarry Smith ierr = MatInterpolate(mg[i+1]->interpolate,mg[i]->x,mg[i+1]->x);CHKERRQ(ierr); 76*4b9ad928SBarry Smith } 77*4b9ad928SBarry Smith if (mg[l-1]->eventsolve) {ierr = PetscLogEventBegin(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);} 78*4b9ad928SBarry Smith ierr = KSPSetRhs(mg[l-1]->smoothd,mg[l-1]->b);CHKERRQ(ierr); 79*4b9ad928SBarry Smith ierr = KSPSetSolution(mg[l-1]->smoothd,mg[l-1]->x);CHKERRQ(ierr); 80*4b9ad928SBarry Smith ierr = KSPSolve(mg[l-1]->smoothd);CHKERRQ(ierr); 81*4b9ad928SBarry Smith if (mg[l-1]->eventsolve) {ierr = PetscLogEventEnd(mg[l-1]->eventsolve,0,0,0,0);CHKERRQ(ierr);} 82*4b9ad928SBarry Smith 83*4b9ad928SBarry Smith PetscFunctionReturn(0); 84*4b9ad928SBarry Smith } 85*4b9ad928SBarry Smith 86*4b9ad928SBarry Smith 87