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 { 9f3fbd535SBarry Smith PetscInt i,l = mglevels[0]->levels; 104b9ad928SBarry Smith 114b9ad928SBarry Smith PetscFunctionBegin; 124b9ad928SBarry Smith /* compute RHS on each level */ 134b9ad928SBarry Smith for (i=l-1; i>0; i--) { 145f80ce2aSJacob Faibussowitsch if (mglevels[i]->eventinterprestrict) CHKERRQ(PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0)); 15fcb023d4SJed Brown if (!transpose) { 165f80ce2aSJacob Faibussowitsch if (matapp) CHKERRQ(MatMatRestrict(mglevels[i]->restrct,mglevels[i]->B,&mglevels[i-1]->B)); 175f80ce2aSJacob Faibussowitsch else CHKERRQ(MatRestrict(mglevels[i]->restrct,mglevels[i]->b,mglevels[i-1]->b)); 18fcb023d4SJed Brown } else { 195f80ce2aSJacob Faibussowitsch if (matapp) CHKERRQ(MatMatRestrict(mglevels[i]->interpolate,mglevels[i]->B,&mglevels[i-1]->B)); 205f80ce2aSJacob Faibussowitsch else CHKERRQ(MatRestrict(mglevels[i]->interpolate,mglevels[i]->b,mglevels[i-1]->b)); 21fcb023d4SJed Brown } 225f80ce2aSJacob Faibussowitsch if (mglevels[i]->eventinterprestrict) CHKERRQ(PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0)); 234b9ad928SBarry Smith } 24a8c7a070SBarry Smith /* solve separately on each level */ 254b9ad928SBarry Smith for (i=0; i<l; i++) { 2630b0564aSStefano Zampini if (matapp) { 2730b0564aSStefano Zampini if (!mglevels[i]->X) { 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(mglevels[i]->B,MAT_DO_NOT_COPY_VALUES,&mglevels[i]->X)); 2930b0564aSStefano Zampini } else { 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(mglevels[i]->X)); 3130b0564aSStefano Zampini } 3230b0564aSStefano Zampini } else { 335f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(mglevels[i]->x)); 3430b0564aSStefano Zampini } 355f80ce2aSJacob Faibussowitsch if (mglevels[i]->eventsmoothsolve) CHKERRQ(PetscLogEventBegin(mglevels[i]->eventsmoothsolve,0,0,0,0)); 36fcb023d4SJed Brown if (!transpose) { 3730b0564aSStefano Zampini if (matapp) { 385f80ce2aSJacob Faibussowitsch CHKERRQ(KSPMatSolve(mglevels[i]->smoothd,mglevels[i]->B,mglevels[i]->X)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCheckSolve(mglevels[i]->smoothd,pc,NULL)); 4030b0564aSStefano Zampini } else { 415f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(mglevels[i]->smoothd,mglevels[i]->b,mglevels[i]->x)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCheckSolve(mglevels[i]->smoothd,pc,mglevels[i]->x)); 4330b0564aSStefano Zampini } 44fcb023d4SJed Brown } else { 45*28b400f6SJacob Faibussowitsch PetscCheck(!matapp,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not supported"); 465f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolveTranspose(mglevels[i]->smoothu,mglevels[i]->b,mglevels[i]->x)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCheckSolve(mglevels[i]->smoothu,pc,mglevels[i]->x)); 48fcb023d4SJed Brown } 495f80ce2aSJacob Faibussowitsch if (mglevels[i]->eventsmoothsolve) CHKERRQ(PetscLogEventEnd(mglevels[i]->eventsmoothsolve,0,0,0,0)); 504b9ad928SBarry Smith } 514b9ad928SBarry Smith for (i=1; i<l; i++) { 525f80ce2aSJacob Faibussowitsch if (mglevels[i]->eventinterprestrict) CHKERRQ(PetscLogEventBegin(mglevels[i]->eventinterprestrict,0,0,0,0)); 53fcb023d4SJed Brown if (!transpose) { 545f80ce2aSJacob Faibussowitsch if (matapp) CHKERRQ(MatMatInterpolateAdd(mglevels[i]->interpolate,mglevels[i-1]->X,mglevels[i]->X,&mglevels[i]->X)); 555f80ce2aSJacob Faibussowitsch else CHKERRQ(MatInterpolateAdd(mglevels[i]->interpolate,mglevels[i-1]->x,mglevels[i]->x,mglevels[i]->x)); 56fcb023d4SJed Brown } else { 575f80ce2aSJacob Faibussowitsch if (matapp) CHKERRQ(MatMatInterpolateAdd(mglevels[i]->restrct,mglevels[i-1]->X,mglevels[i]->X,&mglevels[i]->X)); 585f80ce2aSJacob Faibussowitsch else CHKERRQ(MatInterpolateAdd(mglevels[i]->restrct,mglevels[i-1]->x,mglevels[i]->x,mglevels[i]->x)); 59fcb023d4SJed Brown } 605f80ce2aSJacob Faibussowitsch if (mglevels[i]->eventinterprestrict) CHKERRQ(PetscLogEventEnd(mglevels[i]->eventinterprestrict,0,0,0,0)); 614b9ad928SBarry Smith } 624b9ad928SBarry Smith PetscFunctionReturn(0); 634b9ad928SBarry Smith } 64