From 5693bf9c952f893490a745af18cabb8ce32d09ce Mon Sep 17 00:00:00 2001 From: Weiqun Zhang Date: Wed, 20 Dec 2023 17:08:00 -0800 Subject: [PATCH] wip --- Src/LinearSolvers/MLMG/AMReX_MLCurlCurl.H | 42 ++++++++++++++++------- 1 file changed, 29 insertions(+), 13 deletions(-) diff --git a/Src/LinearSolvers/MLMG/AMReX_MLCurlCurl.H b/Src/LinearSolvers/MLMG/AMReX_MLCurlCurl.H index b836f40494..62faf1d47c 100644 --- a/Src/LinearSolvers/MLMG/AMReX_MLCurlCurl.H +++ b/Src/LinearSolvers/MLMG/AMReX_MLCurlCurl.H @@ -120,11 +120,15 @@ public: [[nodiscard]] MF makeCoarseAmr (int famrlev, IntVect const& ng) const override; private: + + void applyBC (int amrlev, int mglev, MF& in, BCMode bc_mode, StateMode s_mode) const; + RT m_alpha = std::numeric_limits::lowest(); RT m_beta = std::numeric_limits::lowest(); Array m_etype{IntVect(0,1,1), IntVect(1,0,1), IntVect(1,1,0)}; + static constexpr int m_ncomp = 1; using mf_type = typename MF::value_type; }; @@ -208,10 +212,11 @@ MLCurlCurlT::averageDownSolutionRHS (int camrlev, MF& crse_sol, MF& crse_rhs template void MLCurlCurlT::apply (int amrlev, int mglev, MF& out, MF& in, BCMode bc_mode, - StateMode s_mode, const MLMGBndryT* bndry) const + StateMode s_mode, const MLMGBndryT* /*bndry*/) const { - amrex::ignore_unused(amrlev, mglev, out, in, bc_mode, s_mode, bndry); - amrex::Abort("MLCurlCurlT::apply: TODO"); + applyBC(amrlev, mglev, in, bc_mode, s_mode); + + amrex::Warning("MLCurlCurlT::apply: TODO"); } template @@ -226,10 +231,14 @@ MLCurlCurlT::smooth (int amrlev, int mglev, MF& sol, const MF& rhs, template void MLCurlCurlT::solutionResidual (int amrlev, MF& resid, MF& x, const MF& b, - const MF* crse_bcdata) + const MF* /*crse_bcdata*/) { - amrex::ignore_unused(amrlev, resid, x, b, crse_bcdata); - amrex::Abort("MLCurlCurlT::solutionResidual: TODO"); + BL_PROFILE("MLCurlCurl::solutionResidual()"); + const int mglev = 0; + apply(amrlev, mglev, resid, x, BCMode::Inhomogeneous, StateMode::Solution); + for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) { + amrex::Xpay(resid, RT(-1.0), b, 0, 0, m_ncomp, IntVect(0)); + } } template @@ -288,11 +297,9 @@ MLCurlCurlT::xdoty (int amrlev, int mglev, const MF& x, const MF& y, bool lo template typename MLCurlCurlT::RT -MLCurlCurlT::normInf (int amrlev, MF const& mf, bool local) const +MLCurlCurlT::normInf (int /*amrlev*/, MF const& mf, bool local) const { - amrex::ignore_unused(amrlev, mf, local); - amrex::Abort("MLCurlCurlT::normInf: TODO"); - return 0; + return amrex::norminf(mf, 0, m_ncomp, IntVect(0), local); } template @@ -333,7 +340,7 @@ MLCurlCurlT::make (int amrlev, int mglev, IntVect const& ng) const MF r; for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) { r[idim].define(amrex::convert(this->m_grids[amrlev][mglev], m_etype[idim]), - this->m_dmap[amrlev][mglev], this->getNComp(), ng, MFInfo(), + this->m_dmap[amrlev][mglev], m_ncomp, ng, MFInfo(), *(this->m_factory)[amrlev][mglev]); } return r; @@ -361,7 +368,7 @@ MLCurlCurlT::makeCoarseMG (int amrlev, int mglev, IntVect const& ng) const MF r; for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) { r[idim].define(amrex::convert(cba, m_etype[idim]), - this->m_dmap[amrlev][mglev], this->getNComp(), ng); + this->m_dmap[amrlev][mglev], m_ncomp, ng); } return r; } @@ -377,11 +384,20 @@ MLCurlCurlT::makeCoarseAmr (int famrlev, IntVect const& ng) const MF r; for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) { r[idim].define(amrex::convert(cba, m_etype[idim]), - this->m_dmap[famrlev][0], this->getNComp(), ng); + this->m_dmap[famrlev][0], m_ncomp, ng); } return r; } +template +void +MLCurlCurlT::applyBC (int amrlev, int mglev, MF& in, BCMode /*bc_mode*/, + StateMode /*s_mode*/) const +{ + Vector mfs{AMREX_D_DECL(&(in[0]),&(in[1]),&(in[2]))}; + FillBoundary(mfs, this->m_geom[amrlev][mglev].periodicity()); +} + template <> struct LinOpData , void> {