14b9ad928SBarry Smith /* 24b9ad928SBarry Smith Full multigrid using either additive or multiplicative V or W cycle 34b9ad928SBarry Smith */ 4af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> 54b9ad928SBarry Smith 6*30b0564aSStefano Zampini PetscErrorCode PCMGFCycle_Private(PC pc,PC_MG_Levels **mglevels,PetscBool transpose,PetscBool matapp) 74b9ad928SBarry Smith { 86849ba73SBarry Smith PetscErrorCode ierr; 9f3fbd535SBarry Smith PetscInt i,l = mglevels[0]->levels; 104b9ad928SBarry Smith 114b9ad928SBarry Smith PetscFunctionBegin; 12fcb023d4SJed Brown if (!transpose) { 134b9ad928SBarry Smith /* restrict the RHS through all levels to coarsest. */ 144b9ad928SBarry Smith for (i=l-1; i>0; i--) { 1563e6d426SJed Brown if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 16*30b0564aSStefano Zampini if (matapp) { ierr = MatMatRestrict(mglevels[i]->restrct,mglevels[i]->B,&mglevels[i-1]->B);CHKERRQ(ierr); } 17*30b0564aSStefano Zampini else { ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr); } 1863e6d426SJed Brown if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 194b9ad928SBarry Smith } 204b9ad928SBarry Smith 214b9ad928SBarry Smith /* work our way up through the levels */ 22*30b0564aSStefano Zampini if (matapp) { 23*30b0564aSStefano Zampini if (!mglevels[0]->X) { 24*30b0564aSStefano Zampini ierr = MatDuplicate(mglevels[0]->B,MAT_DO_NOT_COPY_VALUES,&mglevels[0]->X);CHKERRQ(ierr); 25*30b0564aSStefano Zampini } else { 26*30b0564aSStefano Zampini ierr = MatZeroEntries(mglevels[0]->X);CHKERRQ(ierr); 27*30b0564aSStefano Zampini } 28*30b0564aSStefano Zampini } else { 29*30b0564aSStefano Zampini ierr = VecZeroEntries(mglevels[0]->x);CHKERRQ(ierr); 30*30b0564aSStefano Zampini } 314b9ad928SBarry Smith for (i=0; i<l-1; i++) { 32*30b0564aSStefano Zampini ierr = PCMGMCycle_Private(pc,&mglevels[i],transpose,matapp,NULL);CHKERRQ(ierr); 3363e6d426SJed Brown if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 34*30b0564aSStefano Zampini if (matapp) { ierr = MatMatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->X,&mglevels[i+1]->X);CHKERRQ(ierr); } 35*30b0564aSStefano Zampini else { ierr = MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);CHKERRQ(ierr); } 3663e6d426SJed Brown if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 374b9ad928SBarry Smith } 38*30b0564aSStefano Zampini ierr = PCMGMCycle_Private(pc,&mglevels[l-1],transpose,matapp,NULL);CHKERRQ(ierr); 39fcb023d4SJed Brown } else { 40*30b0564aSStefano Zampini ierr = PCMGMCycle_Private(pc,&mglevels[l-1],transpose,matapp,NULL);CHKERRQ(ierr); 41fcb023d4SJed Brown for (i=l-2; i>=0; i--) { 42fcb023d4SJed Brown if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 43*30b0564aSStefano Zampini if (matapp) { ierr = MatMatRestrict(mglevels[i+1]->interpolate,mglevels[i+1]->X,&mglevels[i]->X);CHKERRQ(ierr); } 44*30b0564aSStefano Zampini else { ierr = MatRestrict(mglevels[i+1]->interpolate,mglevels[i+1]->x,mglevels[i]->x);CHKERRQ(ierr); } 45fcb023d4SJed Brown if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 46*30b0564aSStefano Zampini ierr = PCMGMCycle_Private(pc,&mglevels[i],transpose,matapp,NULL);CHKERRQ(ierr); 47fcb023d4SJed Brown } 48fcb023d4SJed Brown for (i=1; i<l; i++) { 49fcb023d4SJed Brown if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 50*30b0564aSStefano Zampini if (matapp) { ierr = MatMatInterpolate(mglevels[i]->restrct,mglevels[i-1]->B,&mglevels[i]->B);CHKERRQ(ierr); } 51*30b0564aSStefano Zampini else { ierr = MatInterpolate(mglevels[i]->restrct,mglevels[i-1]->b,mglevels[i]->b);CHKERRQ(ierr); } 52fcb023d4SJed Brown if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 53fcb023d4SJed Brown } 54fcb023d4SJed Brown } 554b9ad928SBarry Smith PetscFunctionReturn(0); 564b9ad928SBarry Smith } 574b9ad928SBarry Smith 58*30b0564aSStefano Zampini PetscErrorCode PCMGKCycle_Private(PC pc,PC_MG_Levels **mglevels,PetscBool transpose,PetscBool matapp) 594b9ad928SBarry Smith { 606849ba73SBarry Smith PetscErrorCode ierr; 61f3fbd535SBarry Smith PetscInt i,l = mglevels[0]->levels; 624b9ad928SBarry Smith 634b9ad928SBarry Smith PetscFunctionBegin; 644b9ad928SBarry Smith /* restrict the RHS through all levels to coarsest. */ 654b9ad928SBarry Smith for (i=l-1; i>0; i--) { 6663e6d426SJed Brown if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 67*30b0564aSStefano Zampini if (matapp) { ierr = MatMatRestrict(mglevels[i]->restrct,mglevels[i]->B,&mglevels[i-1]->B);CHKERRQ(ierr); } 68*30b0564aSStefano Zampini else { ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr); } 6963e6d426SJed Brown if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 704b9ad928SBarry Smith } 714b9ad928SBarry Smith 724b9ad928SBarry Smith /* work our way up through the levels */ 73*30b0564aSStefano Zampini if (matapp) { 74*30b0564aSStefano Zampini if (!mglevels[0]->X) { 75*30b0564aSStefano Zampini ierr = MatDuplicate(mglevels[0]->B,MAT_DO_NOT_COPY_VALUES,&mglevels[0]->X);CHKERRQ(ierr); 76*30b0564aSStefano Zampini } else { 77*30b0564aSStefano Zampini ierr = MatZeroEntries(mglevels[0]->X);CHKERRQ(ierr); 78*30b0564aSStefano Zampini } 79*30b0564aSStefano Zampini } else { 80*30b0564aSStefano Zampini ierr = VecZeroEntries(mglevels[0]->x);CHKERRQ(ierr); 81*30b0564aSStefano Zampini } 824b9ad928SBarry Smith for (i=0; i<l-1; i++) { 8363e6d426SJed Brown if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 84*30b0564aSStefano Zampini if (matapp) { 85*30b0564aSStefano Zampini ierr = KSPMatSolve(mglevels[i]->smoothd,mglevels[i]->B,mglevels[i]->X);CHKERRQ(ierr); 86*30b0564aSStefano Zampini ierr = KSPCheckSolve(mglevels[i]->smoothd,pc,NULL);CHKERRQ(ierr); 87*30b0564aSStefano Zampini } else { 88f3fbd535SBarry Smith ierr = KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);CHKERRQ(ierr); 89c0decd05SBarry Smith ierr = KSPCheckSolve(mglevels[i]->smoothd,pc,mglevels[i]->x);CHKERRQ(ierr); 90*30b0564aSStefano Zampini } 9163e6d426SJed Brown if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 9263e6d426SJed Brown if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 93*30b0564aSStefano Zampini if (matapp) { ierr = MatMatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->X,&mglevels[i+1]->X);CHKERRQ(ierr); } 94*30b0564aSStefano Zampini else { ierr = MatInterpolate(mglevels[i+1]->interpolate,mglevels[i]->x,mglevels[i+1]->x);CHKERRQ(ierr); } 9563e6d426SJed Brown if (mglevels[i+1]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i+1]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 964b9ad928SBarry Smith } 9763e6d426SJed Brown if (mglevels[l-1]->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 98*30b0564aSStefano Zampini if (matapp) { 99*30b0564aSStefano Zampini ierr = KSPMatSolve(mglevels[l-1]->smoothd,mglevels[l-1]->B,mglevels[l-1]->X);CHKERRQ(ierr); 100*30b0564aSStefano Zampini ierr = KSPCheckSolve(mglevels[l-1]->smoothd,pc,NULL);CHKERRQ(ierr); 101*30b0564aSStefano Zampini } else { 102f3fbd535SBarry Smith ierr = KSPSolve(mglevels[l-1]->smoothd,mglevels[l-1]->b,mglevels[l-1]->x);CHKERRQ(ierr); 103c0decd05SBarry Smith ierr = KSPCheckSolve(mglevels[l-1]->smoothd,pc,mglevels[l-1]->x);CHKERRQ(ierr); 104*30b0564aSStefano Zampini } 10563e6d426SJed Brown if (mglevels[l-1]->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels[l-1]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 1064b9ad928SBarry Smith PetscFunctionReturn(0); 1074b9ad928SBarry Smith } 108