1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith Additive Multigrid V Cycle routine 44b9ad928SBarry Smith */ 5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> 64b9ad928SBarry Smith 730b0564aSStefano Zampini PetscErrorCode PCMGACycle_Private(PC pc,PC_MG_Levels **mglevels,PetscBool transpose,PetscBool matapp) 84b9ad928SBarry Smith { 96849ba73SBarry Smith PetscErrorCode ierr; 10f3fbd535SBarry Smith PetscInt i,l = mglevels[0]->levels; 114b9ad928SBarry Smith 124b9ad928SBarry Smith PetscFunctionBegin; 134b9ad928SBarry Smith /* compute RHS on each level */ 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);} 16fcb023d4SJed Brown if (!transpose) { 1730b0564aSStefano Zampini if (matapp) { ierr = MatMatRestrict(mglevels[i]->restrct,mglevels[i]->B,&mglevels[i-1]->B);CHKERRQ(ierr); } 1830b0564aSStefano Zampini else { ierr = MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr); } 19fcb023d4SJed Brown } else { 2030b0564aSStefano Zampini if (matapp) { ierr = MatMatRestrict(mglevels[i]->interpolate,mglevels[i]->B,&mglevels[i-1]->B);CHKERRQ(ierr); } 2130b0564aSStefano Zampini else { ierr = MatRestrict(mglevels[i]->interpolate,mglevels[i]->b,mglevels[i-1]->b);CHKERRQ(ierr); } 22fcb023d4SJed Brown } 2363e6d426SJed Brown if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 244b9ad928SBarry Smith } 25a8c7a070SBarry Smith /* solve separately on each level */ 264b9ad928SBarry Smith for (i=0; i<l; i++) { 2730b0564aSStefano Zampini if (matapp) { 2830b0564aSStefano Zampini if (!mglevels[i]->X) { 2930b0564aSStefano Zampini ierr = MatDuplicate(mglevels[i]->B,MAT_DO_NOT_COPY_VALUES,&mglevels[i]->X);CHKERRQ(ierr); 3030b0564aSStefano Zampini } else { 3130b0564aSStefano Zampini ierr = MatZeroEntries(mglevels[i]->X);CHKERRQ(ierr); 3230b0564aSStefano Zampini } 3330b0564aSStefano Zampini } else { 3430b0564aSStefano Zampini ierr = VecZeroEntries(mglevels[i]->x);CHKERRQ(ierr); 3530b0564aSStefano Zampini } 3663e6d426SJed Brown if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 37fcb023d4SJed Brown if (!transpose) { 3830b0564aSStefano Zampini if (matapp) { 3930b0564aSStefano Zampini ierr = KSPMatSolve(mglevels[i]->smoothd,mglevels[i]->B,mglevels[i]->X);CHKERRQ(ierr); 4030b0564aSStefano Zampini ierr = KSPCheckSolve(mglevels[i]->smoothd,pc,NULL);CHKERRQ(ierr); 4130b0564aSStefano Zampini } else { 42f3fbd535SBarry Smith ierr = KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x);CHKERRQ(ierr); 43c0decd05SBarry Smith ierr = KSPCheckSolve(mglevels[i]->smoothd,pc,mglevels[i]->x);CHKERRQ(ierr); 4430b0564aSStefano Zampini } 45fcb023d4SJed Brown } else { 46*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 47fcb023d4SJed Brown ierr = KSPSolveTranspose(mglevels[i]->smoothu,mglevels[i]->b,mglevels[i]->x);CHKERRQ(ierr); 48fcb023d4SJed Brown ierr = KSPCheckSolve(mglevels[i]->smoothu,pc,mglevels[i]->x);CHKERRQ(ierr); 49fcb023d4SJed Brown } 5063e6d426SJed Brown if (mglevels[i]->eventsmoothsolve) {ierr = PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 514b9ad928SBarry Smith } 524b9ad928SBarry Smith for (i=1; i<l; i++) { 5363e6d426SJed Brown if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 54fcb023d4SJed Brown if (!transpose) { 5530b0564aSStefano Zampini if (matapp) { ierr = MatMatInterpolateAdd(mglevels[i]->interpolate,mglevels[i-1]->X,mglevels[i]->X,&mglevels[i]->X);CHKERRQ(ierr); } 5630b0564aSStefano Zampini else { ierr = MatInterpolateAdd(mglevels[i]->interpolate,mglevels[i-1]->x,mglevels[i]->x,mglevels[i]->x);CHKERRQ(ierr); } 57fcb023d4SJed Brown } else { 5830b0564aSStefano Zampini if (matapp) { ierr = MatMatInterpolateAdd(mglevels[i]->restrct,mglevels[i-1]->X,mglevels[i]->X,&mglevels[i]->X);CHKERRQ(ierr); } 5930b0564aSStefano Zampini else { ierr = MatInterpolateAdd(mglevels[i]->restrct,mglevels[i-1]->x,mglevels[i]->x,mglevels[i]->x);CHKERRQ(ierr); } 60fcb023d4SJed Brown } 6163e6d426SJed Brown if (mglevels[i]->eventinterprestrict) {ierr = PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 624b9ad928SBarry Smith } 634b9ad928SBarry Smith PetscFunctionReturn(0); 644b9ad928SBarry Smith } 65