Skip to content

Commit

Permalink
Fix MLMG::getGradSolution & getFluxes for inhomogeneous Neumann and R…
Browse files Browse the repository at this point in the history
…obin BC

Because of the way how inhomogeneous and Robin BC are handled, we must add
the inhomogeneous fluxes back, otherwise they would be zero at those
boundaries.
  • Loading branch information
WeiqunZhang committed Oct 11, 2022
1 parent 5acfe07 commit 76e058b
Show file tree
Hide file tree
Showing 5 changed files with 164 additions and 0 deletions.
4 changes: 4 additions & 0 deletions Src/LinearSolvers/MLMG/AMReX_MLCellABecLap.H
Original file line number Diff line number Diff line change
Expand Up @@ -61,6 +61,10 @@ public:

virtual void applyInhomogNeumannTerm (int amrlev, Any& rhs) const final override;

virtual void addInhomogNeumannFlux (
int amrlev, const Array<MultiFab*,AMREX_SPACEDIM>& grad,
MultiFab const& sol, bool mult_bcoef) const final override;

virtual void applyOverset (int amlev, Any& rhs) const override;

#if defined(AMREX_USE_HYPRE) && (AMREX_SPACEDIM > 1)
Expand Down
111 changes: 111 additions & 0 deletions Src/LinearSolvers/MLMG/AMReX_MLCellABecLap.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -189,6 +189,7 @@ MLCellABecLap::getFluxes (const Vector<Array<MultiFab*,AMREX_SPACEDIM> >& a_flux
a_flux[alev][idim]->mult(betainv);
}
}
addInhomogNeumannFlux(alev, a_flux[alev], *a_sol[alev], true);
}
}

Expand Down Expand Up @@ -416,6 +417,116 @@ MLCellABecLap::applyInhomogNeumannTerm (int amrlev, Any& a_rhs) const
}
}

void
MLCellABecLap::addInhomogNeumannFlux (
int amrlev, const Array<MultiFab*,AMREX_SPACEDIM>& grad, MultiFab const& sol,
bool mult_bcoef) const
{
/*
* if (mult_bcoef == true)
* grad is -bceof*grad phi
* else
* grad is grad phi
*/
Real fac = mult_bcoef ? Real(-1.0) : Real(1.0);

bool has_inhomog_neumann = hasInhomogNeumannBC();
bool has_robin = hasRobinBC();

if (!has_inhomog_neumann && !has_robin) return;

int ncomp = getNComp();
const int mglev = 0;

const auto dxinv = m_geom[amrlev][mglev].InvCellSize();
const Box domain = m_geom[amrlev][mglev].growPeriodicDomain(1);

Array<MultiFab const*, AMREX_SPACEDIM> bcoef = {AMREX_D_DECL(nullptr,nullptr,nullptr)};
if (mult_bcoef) {
bcoef = getBCoeffs(amrlev,mglev);
}

const auto& bndry = *m_bndry_sol[amrlev];

MFItInfo mfi_info;
if (Gpu::notInLaunchRegion()) mfi_info.SetDynamic(true);

#ifdef AMREX_USE_OMP
#pragma omp parallel if (Gpu::notInLaunchRegion())
#endif
for (MFIter mfi(sol, mfi_info); mfi.isValid(); ++mfi)
{
Box const& vbx = mfi.validbox();
for (OrientationIter orit; orit.isValid(); ++orit) {
const Orientation ori = orit();
const int idim = ori.coordDir();
const Box& ccb = amrex::adjCell(vbx, ori);
const Dim3 os = IntVect::TheDimensionVector(idim).dim3();
const Real dxi = dxinv[idim];
if (! domain.contains(ccb)) {
for (int icomp = 0; icomp < ncomp; ++icomp) {
auto const& phi = sol.const_array(mfi,icomp);
auto const bv = bndry.bndryValues(ori).multiFab().const_array(mfi,icomp);
auto const bc = bcoef[idim] ? bcoef[idim]->const_array(mfi,icomp)
: Array4<Real const>{};
auto const& f = grad[idim]->array(mfi,icomp);
if (ori.isLow()) {
if (m_lobc_orig[icomp][idim] ==
LinOpBCType::inhomogNeumann) {
AMREX_HOST_DEVICE_FOR_3D(ccb, i, j, k,
{
int ii = i+os.x;
int jj = j+os.y;
int kk = k+os.z;
Real b = bc ? bc(ii,jj,kk) : Real(1.0);
f(ii,jj,kk) = fac*b*bv(i,j,k);
});
} else if (m_lobc_orig[icomp][idim] ==
LinOpBCType::Robin) {
Array4<Real const> const& rbc = (*m_robin_bcval[amrlev])[mfi].const_array(icomp*3);
AMREX_HOST_DEVICE_FOR_3D(ccb, i, j, k,
{
int ii = i+os.x;
int jj = j+os.y;
int kk = k+os.z;
Real tmp = Real(1.0) /
(rbc(i,j,k,1)*dxi + rbc(i,j,k,0)*Real(0.5));
Real RA = rbc(i,j,k,2) * tmp;
Real RB = (rbc(i,j,k,1)*dxi - rbc(i,j,k,0)*Real(0.5)) * tmp;
Real b = bc ? bc(ii,jj,kk) : Real(1.0);
f(ii,jj,kk) = fac*b*dxi*((Real(1.0)-RB)*phi(ii,jj,kk)-RA);
});
}
} else {
if (m_hibc_orig[icomp][idim] ==
LinOpBCType::inhomogNeumann) {
AMREX_HOST_DEVICE_FOR_3D(ccb, i, j, k,
{
Real b = bc ? bc(i,j,k) : Real(1.0);
f(i,j,k) = fac*b*bv(i,j,k);
});
} else if (m_hibc_orig[icomp][idim] ==
LinOpBCType::Robin) {
Array4<Real const> const& rbc = (*m_robin_bcval[amrlev])[mfi].const_array(icomp*3);
AMREX_HOST_DEVICE_FOR_3D(ccb, i, j, k,
{
Real tmp = Real(1.0) /
(rbc(i,j,k,1)*dxi + rbc(i,j,k,0)*Real(0.5));
Real RA = rbc(i,j,k,2) * tmp;
Real RB = (rbc(i,j,k,1)*dxi - rbc(i,j,k,0)*Real(0.5)) * tmp;
Real b = bc ? bc(i,j,k) : Real(1.0);
f(i,j,k) = fac*b*dxi*(RA+(RB-Real(1.0))*
phi(i-os.x,j-os.y,k-os.z));
});
}
}
}
}
}
}
}


void
MLCellABecLap::applyOverset (int amrlev, Any& a_rhs) const
{
Expand Down
42 changes: 42 additions & 0 deletions Src/LinearSolvers/MLMG/AMReX_MLCellABecLap_K.H
Original file line number Diff line number Diff line change
Expand Up @@ -12,4 +12,46 @@
#include <AMReX_MLCellABecLap_3D_K.H>
#endif

namespace amrex {

#if 0
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void mlcellabec_add_innu_xlo (int i, int j, int k, int n, Array<Real> const& fx,
Array<Real const> const& bcoef,
Array4<Real const> const& bcval) noexcept
{
Real b = bcoef ? bcoef(i,j,k,n) : Real(1.0);
fx(i,j,k,n) = -b*bcval(i-1,j,k,icomp);
}

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void mlcellabec_add_innu_ylo (int i, int j, int k, int n, Array<Real> const& fy,
Array<Real const> const& bcoef,
Array4<Real const> const& bcval) noexcept
{
Real b = bcoef ? bcoef(i,j,k,n) : Real(1.0);
fy(i,j,k,n) = -b*bcval(i,j-1,k,icomp);
}

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void mlcellabec_add_innu_zlo (int i, int j, int k, int n, Array<Real> const& fz,
Array<Real const> const& bcoef,
Array4<Real const> const& bcval) noexcept
{
Real b = bcoef ? bcoef(i,j,k,n) : Real(1.0);
fz(i,j,k,n) = -b*bcval(i,j,k-1,icomp);
}

AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void mlcellabec_add_innu_hi (int i, int j, int k, int n, Array<Real> const& f,
Array<Real const> const& bcoef,
Array4<Real const> const& bcval) noexcept
{
Real b = bcoef ? bcoef(i,j,k,n) : Real(1.0);
f(i,j,k,n) = -b*bcval(i,j,k,icomp);
}
#endif

}

#endif
5 changes: 5 additions & 0 deletions Src/LinearSolvers/MLMG/AMReX_MLCellLinOp.H
Original file line number Diff line number Diff line change
Expand Up @@ -164,6 +164,11 @@ public:

virtual void AnyAverageDownAndSync (Vector<Any>& sol) const override;

virtual void addInhomogNeumannFlux (int /*amrlev*/,
const Array<MultiFab*,AMREX_SPACEDIM>& /*grad*/,
MultiFab const& /*sol*/,
bool /*mult_bcoef*/) const {}

struct BCTL {
BoundCond type;
Real location;
Expand Down
2 changes: 2 additions & 0 deletions Src/LinearSolvers/MLMG/AMReX_MLCellLinOp.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -938,6 +938,8 @@ MLCellLinOp::compGrad (int amrlev, const Array<MultiFab*,AMREX_SPACEDIM>& grad,
});
#endif
}

addInhomogNeumannFlux(amrlev, grad, sol, false);
}

void
Expand Down

0 comments on commit 76e058b

Please sign in to comment.