xref: /petsc/src/ksp/pc/impls/mg/smg.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
34b9ad928SBarry Smith      Additive Multigrid V Cycle routine
44b9ad928SBarry Smith */
5af0996ceSBarry Smith #include <petsc/private/pcmgimpl.h>
64b9ad928SBarry Smith 
7*d71ae5a4SJacob Faibussowitsch PetscErrorCode PCMGACycle_Private(PC pc, PC_MG_Levels **mglevels, PetscBool transpose, PetscBool matapp)
8*d71ae5a4SJacob Faibussowitsch {
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--) {
149566063dSJacob Faibussowitsch     if (mglevels[i]->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels[i]->eventinterprestrict, 0, 0, 0, 0));
15fcb023d4SJed Brown     if (!transpose) {
169566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatRestrict(mglevels[i]->restrct, mglevels[i]->B, &mglevels[i - 1]->B));
179566063dSJacob Faibussowitsch       else PetscCall(MatRestrict(mglevels[i]->restrct, mglevels[i]->b, mglevels[i - 1]->b));
18fcb023d4SJed Brown     } else {
199566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatRestrict(mglevels[i]->interpolate, mglevels[i]->B, &mglevels[i - 1]->B));
209566063dSJacob Faibussowitsch       else PetscCall(MatRestrict(mglevels[i]->interpolate, mglevels[i]->b, mglevels[i - 1]->b));
21fcb023d4SJed Brown     }
229566063dSJacob Faibussowitsch     if (mglevels[i]->eventinterprestrict) PetscCall(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) {
289566063dSJacob Faibussowitsch         PetscCall(MatDuplicate(mglevels[i]->B, MAT_DO_NOT_COPY_VALUES, &mglevels[i]->X));
2930b0564aSStefano Zampini       } else {
309566063dSJacob Faibussowitsch         PetscCall(MatZeroEntries(mglevels[i]->X));
3130b0564aSStefano Zampini       }
3230b0564aSStefano Zampini     } else {
339566063dSJacob Faibussowitsch       PetscCall(VecZeroEntries(mglevels[i]->x));
3430b0564aSStefano Zampini     }
359566063dSJacob Faibussowitsch     if (mglevels[i]->eventsmoothsolve) PetscCall(PetscLogEventBegin(mglevels[i]->eventsmoothsolve, 0, 0, 0, 0));
36fcb023d4SJed Brown     if (!transpose) {
3730b0564aSStefano Zampini       if (matapp) {
389566063dSJacob Faibussowitsch         PetscCall(KSPMatSolve(mglevels[i]->smoothd, mglevels[i]->B, mglevels[i]->X));
399566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(mglevels[i]->smoothd, pc, NULL));
4030b0564aSStefano Zampini       } else {
419566063dSJacob Faibussowitsch         PetscCall(KSPSolve(mglevels[i]->smoothd, mglevels[i]->b, mglevels[i]->x));
429566063dSJacob Faibussowitsch         PetscCall(KSPCheckSolve(mglevels[i]->smoothd, pc, mglevels[i]->x));
4330b0564aSStefano Zampini       }
44fcb023d4SJed Brown     } else {
4528b400f6SJacob Faibussowitsch       PetscCheck(!matapp, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not supported");
469566063dSJacob Faibussowitsch       PetscCall(KSPSolveTranspose(mglevels[i]->smoothu, mglevels[i]->b, mglevels[i]->x));
479566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(mglevels[i]->smoothu, pc, mglevels[i]->x));
48fcb023d4SJed Brown     }
499566063dSJacob Faibussowitsch     if (mglevels[i]->eventsmoothsolve) PetscCall(PetscLogEventEnd(mglevels[i]->eventsmoothsolve, 0, 0, 0, 0));
504b9ad928SBarry Smith   }
514b9ad928SBarry Smith   for (i = 1; i < l; i++) {
529566063dSJacob Faibussowitsch     if (mglevels[i]->eventinterprestrict) PetscCall(PetscLogEventBegin(mglevels[i]->eventinterprestrict, 0, 0, 0, 0));
53fcb023d4SJed Brown     if (!transpose) {
549566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatInterpolateAdd(mglevels[i]->interpolate, mglevels[i - 1]->X, mglevels[i]->X, &mglevels[i]->X));
559566063dSJacob Faibussowitsch       else PetscCall(MatInterpolateAdd(mglevels[i]->interpolate, mglevels[i - 1]->x, mglevels[i]->x, mglevels[i]->x));
56fcb023d4SJed Brown     } else {
579566063dSJacob Faibussowitsch       if (matapp) PetscCall(MatMatInterpolateAdd(mglevels[i]->restrct, mglevels[i - 1]->X, mglevels[i]->X, &mglevels[i]->X));
589566063dSJacob Faibussowitsch       else PetscCall(MatInterpolateAdd(mglevels[i]->restrct, mglevels[i - 1]->x, mglevels[i]->x, mglevels[i]->x));
59fcb023d4SJed Brown     }
609566063dSJacob Faibussowitsch     if (mglevels[i]->eventinterprestrict) PetscCall(PetscLogEventEnd(mglevels[i]->eventinterprestrict, 0, 0, 0, 0));
614b9ad928SBarry Smith   }
624b9ad928SBarry Smith   PetscFunctionReturn(0);
634b9ad928SBarry Smith }
64