14b9ad928SBarry Smith /* 24b9ad928SBarry Smith Additive Multigrid V Cycle routine 34b9ad928SBarry Smith */ 4af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h> 54b9ad928SBarry Smith 6d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGACycle_Private(PC pc, PC_MG_Levels **mglevels, PetscBool transpose, PetscBool matapp) 7d71ae5a4SJacob Faibussowitsch { 8f3fbd535SBarry Smith PetscInt i, l = mglevels[0]->levels; 94b9ad928SBarry Smith 104b9ad928SBarry Smith PetscFunctionBegin; 114b9ad928SBarry Smith /* compute RHS on each level */ 124b9ad928SBarry Smith for (i = l - 1; i > 0; i--) { 139566063dSJacob Faibussowitsch if (mglevels[i]->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels[i]->eventinterprestrict, 0, 0, 0, 0)); 14fcb023d4SJed Brown if (!transpose) { 159566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatRestrict(mglevels[i]->restrct, mglevels[i]->B, &mglevels[i - 1]->B)); 169566063dSJacob Faibussowitsch else PetscCall(MatRestrict(mglevels[i]->restrct, mglevels[i]->b, mglevels[i - 1]->b)); 17fcb023d4SJed Brown } else { 189566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatRestrict(mglevels[i]->interpolate, mglevels[i]->B, &mglevels[i - 1]->B)); 199566063dSJacob Faibussowitsch else PetscCall(MatRestrict(mglevels[i]->interpolate, mglevels[i]->b, mglevels[i - 1]->b)); 20fcb023d4SJed Brown } 219566063dSJacob Faibussowitsch if (mglevels[i]->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels[i]->eventinterprestrict, 0, 0, 0, 0)); 224b9ad928SBarry Smith } 23a8c7a070SBarry Smith /* solve separately on each level */ 244b9ad928SBarry Smith for (i = 0; i < l; i++) { 2530b0564aSStefano Zampini if (matapp) { 2630b0564aSStefano Zampini if (!mglevels[i]->X) { 279566063dSJacob Faibussowitsch PetscCall(MatDuplicate(mglevels[i]->B, MAT_DO_NOT_COPY_VALUES, &mglevels[i]->X)); 2830b0564aSStefano Zampini } else { 299566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mglevels[i]->X)); 3030b0564aSStefano Zampini } 3130b0564aSStefano Zampini } else { 329566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(mglevels[i]->x)); 3330b0564aSStefano Zampini } 349566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsolve, 0, 0, 0, 0)); 35fcb023d4SJed Brown if (!transpose) { 3630b0564aSStefano Zampini if (matapp) { 379566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(mglevels[i]->smoothd, mglevels[i]->B, mglevels[i]->X)); 389566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels[i]->smoothd, pc, NULL)); 3930b0564aSStefano Zampini } else { 409566063dSJacob Faibussowitsch PetscCall(KSPSolve(mglevels[i]->smoothd, mglevels[i]->b, mglevels[i]->x)); 419566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels[i]->smoothd, pc, mglevels[i]->x)); 4230b0564aSStefano Zampini } 43fcb023d4SJed Brown } else { 4428b400f6SJacob Faibussowitsch PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported"); 459566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(mglevels[i]->smoothu, mglevels[i]->b, mglevels[i]->x)); 469566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(mglevels[i]->smoothu, pc, mglevels[i]->x)); 47fcb023d4SJed Brown } 489566063dSJacob Faibussowitsch if (mglevels[i]->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsolve, 0, 0, 0, 0)); 494b9ad928SBarry Smith } 504b9ad928SBarry Smith for (i = 1; i < l; i++) { 519566063dSJacob Faibussowitsch if (mglevels[i]->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels[i]->eventinterprestrict, 0, 0, 0, 0)); 52fcb023d4SJed Brown if (!transpose) { 539566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatInterpolateAdd(mglevels[i]->interpolate, mglevels[i - 1]->X, mglevels[i]->X, &mglevels[i]->X)); 549566063dSJacob Faibussowitsch else PetscCall(MatInterpolateAdd(mglevels[i]->interpolate, mglevels[i - 1]->x, mglevels[i]->x, mglevels[i]->x)); 55fcb023d4SJed Brown } else { 569566063dSJacob Faibussowitsch if (matapp) PetscCall(MatMatInterpolateAdd(mglevels[i]->restrct, mglevels[i - 1]->X, mglevels[i]->X, &mglevels[i]->X)); 579566063dSJacob Faibussowitsch else PetscCall(MatInterpolateAdd(mglevels[i]->restrct, mglevels[i - 1]->x, mglevels[i]->x, mglevels[i]->x)); 58fcb023d4SJed Brown } 599566063dSJacob Faibussowitsch if (mglevels[i]->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels[i]->eventinterprestrict, 0, 0, 0, 0)); 604b9ad928SBarry Smith } 61*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 624b9ad928SBarry Smith } 63