Skip to content

Commit

Permalink
CpuBndryFuncFab: Face Data (AMReX-Codes#2545)
Browse files Browse the repository at this point in the history
For completeness, fab_filfc is implemented for CpuBndryFuncFab.  However,
one should use GpuBndryFuncFab that works for both GPU and CPU builds.
  • Loading branch information
WeiqunZhang authored Jan 3, 2022
1 parent a73ce67 commit d36950d
Show file tree
Hide file tree
Showing 11 changed files with 123 additions and 9 deletions.
2 changes: 1 addition & 1 deletion Src/Base/AMReX_FilCC_1D_C.H
Original file line number Diff line number Diff line change
Expand Up @@ -10,7 +10,7 @@ namespace amrex {

struct FilccCell
{
AMREX_GPU_DEVICE AMREX_FORCE_INLINE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void operator() (const IntVect& iv, Array4<Real> const& q,
const int dcomp, const int numcomp,
Box const& domain_box, const BCRec* bcr,
Expand Down
2 changes: 1 addition & 1 deletion Src/Base/AMReX_FilCC_2D_C.H
Original file line number Diff line number Diff line change
Expand Up @@ -10,7 +10,7 @@ namespace amrex {

struct FilccCell
{
AMREX_GPU_DEVICE AMREX_FORCE_INLINE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void operator() (const IntVect& iv, Array4<Real> const& q,
const int dcomp, const int numcomp,
Box const& domain_box, const BCRec* bcr,
Expand Down
2 changes: 1 addition & 1 deletion Src/Base/AMReX_FilCC_3D_C.H
Original file line number Diff line number Diff line change
Expand Up @@ -10,7 +10,7 @@ namespace amrex {

struct FilccCell
{
AMREX_GPU_DEVICE AMREX_FORCE_INLINE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void operator() (const IntVect& iv, Array4<Real> const& q,
const int dcomp, const int numcomp,
Box const& domain_box, const BCRec* bcr,
Expand Down
3 changes: 1 addition & 2 deletions Src/Base/AMReX_FilFC_1D_C.H
Original file line number Diff line number Diff line change
Expand Up @@ -10,7 +10,7 @@ namespace amrex {

struct FilfcFace
{
AMREX_GPU_DEVICE AMREX_FORCE_INLINE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void operator() (const IntVect& iv, Array4<Real> const& q,
const int dcomp, const int numcomp,
Box const& domain_box, const BCRec* bcr,
Expand Down Expand Up @@ -115,4 +115,3 @@ struct FilfcFace
}

#endif

2 changes: 1 addition & 1 deletion Src/Base/AMReX_FilFC_2D_C.H
Original file line number Diff line number Diff line change
Expand Up @@ -10,7 +10,7 @@ namespace amrex {

struct FilfcFace
{
AMREX_GPU_DEVICE AMREX_FORCE_INLINE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void operator() (const IntVect& iv, Array4<Real> const& q,
const int dcomp, const int numcomp,
Box const& domain_box, const BCRec* bcr,
Expand Down
2 changes: 1 addition & 1 deletion Src/Base/AMReX_FilFC_3D_C.H
Original file line number Diff line number Diff line change
Expand Up @@ -10,7 +10,7 @@ namespace amrex {

struct FilfcFace
{
AMREX_GPU_DEVICE AMREX_FORCE_INLINE
AMREX_GPU_HOST_DEVICE AMREX_FORCE_INLINE
void operator() (const IntVect& iv, Array4<Real> const& q,
const int dcomp, const int numcomp,
Box const& domain_box, const BCRec* bcr,
Expand Down
7 changes: 7 additions & 0 deletions Src/Base/AMReX_FilFC_C.H
Original file line number Diff line number Diff line change
Expand Up @@ -9,4 +9,11 @@
#else
#include <AMReX_FilFC_3D_C.H>
#endif

namespace amrex {
void fab_filfc (Box const& bx, Array4<Real> const& q, int ncomp,
Box const& domain, Real const* dx, Real const* xlo,
BCRec const* bc);
}

#endif
106 changes: 106 additions & 0 deletions Src/Base/AMReX_FilFC_C.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,106 @@
#include <AMReX_FilFC_C.H>

namespace amrex {

void fab_filfc (Box const& bx, Array4<Real> const& qn, int ncomp,
Box const& domain, Real const* /*dx*/, Real const* /*xlo*/,
BCRec const* bcn)
{
Box gdomain = domain;
const auto idxType = bx.ixType();
const IntVect& len = bx.length();
for (int idim = 0; idim < AMREX_SPACEDIM; ++idim) {
if (bcn->lo(idim) == BCType::int_dir) {
AMREX_ASSERT(bcn->hi(idim) == BCType::int_dir);
gdomain.grow(idim, len[idim]);
}
}

FilfcFace fillfunc{};

// fill on the box faces first
{
Array<Box,2*AMREX_SPACEDIM> dom_face_boxes
= {{ AMREX_D_DECL(amrex::convert(amrex::adjCellLo(gdomain, 0, len[0]),idxType),
amrex::convert(amrex::adjCellLo(gdomain, 1, len[1]),idxType),
amrex::convert(amrex::adjCellLo(gdomain, 2, len[2]),idxType)),
AMREX_D_DECL(amrex::convert(amrex::adjCellHi(gdomain, 0, len[0]),idxType),
amrex::convert(amrex::adjCellHi(gdomain, 1, len[1]),idxType),
amrex::convert(amrex::adjCellHi(gdomain, 2, len[2]),idxType)) }};

for (const Box& b : dom_face_boxes) {
Box tmp = b & bx;
amrex::LoopOnCpu(tmp, [=] (int i, int j, int k) noexcept
{
amrex::ignore_unused(j,k);
IntVect const idx(AMREX_D_DECL(i,j,k));
fillfunc(idx, qn, 0, ncomp, domain, bcn, 0);
});
}
}

#if (AMREX_SPACEDIM >= 2)
// fill on the box edges
{
#if (AMREX_SPACEDIM == 2)
Array<Box,4> dom_edge_boxes
= {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType) }};
#else
Array<Box,12> dom_edge_boxes
= {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),idxType),
//
amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),2,len[2]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),2,len[2]),idxType),
//
amrex::convert(amrex::adjCellLo(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellLo(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellLo(gdomain,1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellHi(gdomain,1,len[1]),2,len[2]),idxType) }};
#endif

for (const Box& b : dom_edge_boxes) {
Box tmp = b & bx;
amrex::LoopOnCpu(tmp, [=] (int i, int j, int k) noexcept
{
amrex::ignore_unused(j,k);
IntVect const idx(AMREX_D_DECL(i,j,k));
fillfunc(idx, qn, 0, ncomp, domain, bcn, 0);
});
}
}
#endif

#if (AMREX_SPACEDIM == 3)
// fill on box corners
{
Array<Box,8> dom_corner_boxes
= {{ amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellLo(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellLo(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellLo(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellLo(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType),
amrex::convert(amrex::adjCellHi(amrex::adjCellHi(amrex::adjCellHi(gdomain,0,len[0]),1,len[1]),2,len[2]),idxType) }};

for (const Box& b : dom_corner_boxes) {
Box tmp = b & bx;
amrex::LoopOnCpu(tmp, [=] (int i, int j, int k) noexcept
{
IntVect const idx(AMREX_D_DECL(i,j,k));
fillfunc(idx, qn, 0, ncomp, domain, bcn, 0);
});
}
}
#endif
}

}
2 changes: 1 addition & 1 deletion Src/Base/AMReX_PhysBCFunct.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -66,7 +66,7 @@ CpuBndryFuncFab::operator() (Box const& bx, FArrayBox& dest,
} else if (bx.ixType().nodeCentered()) {
fab_filnd(bx, dest.array(dcomp), numcomp, domain, dx, xlo, &(bcr[bcomp]));
} else {
amrex::Abort("CpuBndryFuncFab: mixed index types are not supported");
fab_filfc(bx, dest.array(dcomp), numcomp, domain, dx, xlo, &(bcr[bcomp]));
}

if (f_user != nullptr)
Expand Down
3 changes: 2 additions & 1 deletion Src/Base/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -162,9 +162,10 @@ target_sources( amrex
AMReX_BC_TYPES.H
AMReX_FilCC_${AMReX_SPACEDIM}D_C.H
AMReX_FilCC_C.H
AMReX_FilCC_C.cpp
AMReX_FilFC_${AMReX_SPACEDIM}D_C.H
AMReX_FilFC_C.H
AMReX_FilCC_C.cpp
AMReX_FilFC_C.cpp
AMReX_FilND_C.H
AMReX_FilND_C.cpp
# Non-Local BC
Expand Down
1 change: 1 addition & 0 deletions Src/Base/Make.package
Original file line number Diff line number Diff line change
Expand Up @@ -178,6 +178,7 @@ C$(AMREX_BASE)_sources += AMReX_FilCC_C.cpp
C$(AMREX_BASE)_headers += AMReX_FilCC_C.H AMReX_FilCC_$(DIM)D_C.H
C$(AMREX_BASE)_sources += AMReX_FilND_C.cpp
C$(AMREX_BASE)_headers += AMReX_FilND_C.H
C$(AMREX_BASE)_sources += AMReX_FilFC_C.cpp
C$(AMREX_BASE)_headers += AMReX_FilFC_C.H AMReX_FilFC_$(DIM)D_C.H

#
Expand Down

0 comments on commit d36950d

Please sign in to comment.