Skip to content

Ek flux boundary continuity #5064

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Open
wants to merge 5 commits into
base: python
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
57 changes: 52 additions & 5 deletions src/walberla_bridge/src/electrokinetics/EKinWalberlaImpl.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -30,6 +30,7 @@
#include <stencil/D3Q27.h>

#include "../BoundaryHandling.hpp"
#include "../BoundaryPackInfo.hpp"
#include "../utils/boundary.hpp"
#include "../utils/types_conversion.hpp"
#include "ek_kernels.hpp"
Expand Down Expand Up @@ -90,6 +91,14 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
using BoundaryModelFlux = BoundaryHandling<Vector3<FloatType>, FixedFlux>;

using BlockStorage = LatticeWalberla::Lattice_T;
struct GhostComm {
/** @brief Ghost communication operations. */
enum GhostCommFlags : unsigned {
FLB, ///< flux boundary communication
DENS, ///< density communication
SIZE
};
};

public:
template <typename T> FloatType FloatType_c(T t) {
Expand Down Expand Up @@ -133,7 +142,7 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
std::shared_ptr<LatticeWalberla> m_lattice;

std::unique_ptr<BoundaryModelDensity> m_boundary_density;
std::unique_ptr<BoundaryModelFlux> m_boundary_flux;
std::shared_ptr<BoundaryModelFlux> m_boundary_flux;

std::unique_ptr<DiffusiveFluxKernel> m_diffusive_flux;
std::unique_ptr<DiffusiveFluxKernelElectrostatic>
Expand Down Expand Up @@ -167,13 +176,18 @@ class EKinWalberlaImpl : public EKinWalberlaBase {

void
reset_flux_boundary_handling(std::shared_ptr<BlockStorage> const &blocks) {
m_boundary_flux = std::make_unique<BoundaryModelFlux>(
m_boundary_flux = std::make_shared<BoundaryModelFlux>(
blocks, m_flux_field_id, m_flag_field_flux_id);
}

using FullCommunicator = blockforest::communication::UniformBufferedScheme<
typename stencil::D3Q27>;
using BoundaryFullCommunicator =
blockforest::communication::UniformBufferedScheme<
typename stencil::D3Q27>;
std::shared_ptr<FullCommunicator> m_full_communication;
std::shared_ptr<BoundaryFullCommunicator> m_boundary_communicator;
std::bitset<GhostComm::SIZE> m_pending_ghost_comm;

public:
EKinWalberlaImpl(std::shared_ptr<LatticeWalberla> lattice, double diffusion,
Expand Down Expand Up @@ -222,6 +236,18 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
m_full_communication->addPackInfo(
std::make_shared<field::communication::PackInfo<DensityField>>(
m_density_field_id));
m_boundary_communicator =
std::make_shared<BoundaryFullCommunicator>(blocks);
m_boundary_communicator->addPackInfo(
std::make_shared<field::communication::PackInfo<FlagField>>(
m_flag_field_flux_id));
auto flux_boundary_packinfo = std::make_shared<
field::communication::BoundaryPackInfo<FlagField, BoundaryModelFlux>>(
m_flag_field_flux_id);
flux_boundary_packinfo->setup_boundary_handle(m_lattice, m_boundary_flux);
m_boundary_communicator->addPackInfo(flux_boundary_packinfo);

m_pending_ghost_comm.set();
}

// Global parameters
Expand Down Expand Up @@ -313,7 +339,20 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
*m_diffusive_flux_electrostatic);
}

void ghost_communication() override { (*m_full_communication)(); }
void ghost_communication() override {
if (m_pending_ghost_comm.test(GhostComm::DENS)) {
(*m_full_communication)();
m_pending_ghost_comm.reset(GhostComm::DENS);
}
ghost_communication_boundary();
}

void ghost_communication_boundary() {
if (m_pending_ghost_comm.test(GhostComm::FLB)) {
m_boundary_communicator->communicate();
m_pending_ghost_comm.reset(GhostComm::FLB);
}
}

private:
void set_diffusion_kernels() {
Expand Down Expand Up @@ -502,6 +541,7 @@ class EKinWalberlaImpl : public EKinWalberlaBase {

// is this the expected behavior when reactions are included?
kernel_boundary_density();
m_pending_ghost_comm.set(GhostComm::DENS);

// Handle VTK writers
integrate_vtk_writers();
Expand All @@ -513,6 +553,7 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
}

bool set_node_density(Utils::Vector3i const &node, double density) override {
m_pending_ghost_comm.set(GhostComm::DENS);
auto bc = get_block_and_cell(get_lattice(), node, false);
if (!bc)
return false;
Expand Down Expand Up @@ -566,6 +607,7 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
void set_slice_density(Utils::Vector3i const &lower_corner,
Utils::Vector3i const &upper_corner,
std::vector<double> const &density) override {
m_pending_ghost_comm.set(GhostComm::DENS);
if (auto const ci = get_interval(lower_corner, upper_corner)) {
auto const &lattice = get_lattice();
auto &block = *(lattice.get_blocks()->begin());
Expand All @@ -587,6 +629,7 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
}

void clear_flux_boundaries() override {
m_pending_ghost_comm.set(GhostComm::FLB);
reset_flux_boundary_handling(get_lattice().get_blocks());
}

Expand All @@ -596,19 +639,20 @@ class EKinWalberlaImpl : public EKinWalberlaBase {

bool set_node_flux_boundary(Utils::Vector3i const &node,
Utils::Vector3d const &flux) override {
m_pending_ghost_comm.set(GhostComm::FLB);
auto bc = get_block_and_cell(get_lattice(), node, true);
if (!bc)
return false;

m_boundary_flux->set_node_value_at_boundary(
node, to_vector3<FloatType>(flux), *bc);

return true;
}

[[nodiscard]] std::optional<Utils::Vector3d>
get_node_flux_at_boundary(Utils::Vector3i const &node,
bool consider_ghosts = false) const override {
assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::FLB)));
auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
if (!bc or !m_boundary_flux->node_is_boundary(node))
return std::nullopt;
Expand All @@ -617,12 +661,12 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
}

bool remove_node_from_flux_boundary(Utils::Vector3i const &node) override {
m_pending_ghost_comm.set(GhostComm::FLB);
auto bc = get_block_and_cell(get_lattice(), node, true);
if (!bc)
return false;

m_boundary_flux->remove_node_from_boundary(node, *bc);

return true;
}

Expand Down Expand Up @@ -710,6 +754,7 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
void set_slice_flux_boundary(
Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
std::vector<std::optional<Utils::Vector3d>> const &flux) override {
m_pending_ghost_comm.set(GhostComm::FLB);
if (auto const ci = get_interval(lower_corner, upper_corner)) {
auto const &lattice = get_lattice();
auto const local_offset = std::get<0>(lattice.get_local_grid_range());
Expand Down Expand Up @@ -804,6 +849,7 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
[[nodiscard]] std::optional<bool>
get_node_is_flux_boundary(Utils::Vector3i const &node,
bool consider_ghosts) const override {
assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::FLB)));
auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
if (!bc)
return std::nullopt;
Expand Down Expand Up @@ -835,6 +881,7 @@ class EKinWalberlaImpl : public EKinWalberlaBase {
void update_flux_boundary_from_shape(
const std::vector<int> &raster_flat,
const std::vector<double> &data_flat) override {
m_pending_ghost_comm.set(GhostComm::FLB);
auto const grid_size = get_lattice().get_grid_dimensions();
auto const data = fill_3D_vector_array(data_flat, grid_size);
set_boundary_from_grid(*m_boundary_flux, get_lattice(), raster_flat, data);
Expand Down
6 changes: 6 additions & 0 deletions src/walberla_bridge/tests/EKinWalberlaImpl_unit_tests.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -161,6 +161,7 @@ BOOST_DATA_TEST_CASE(node_flux_boundary, bdata::make(all_eks()), ek_generator) {
}
{
BOOST_CHECK(ek->set_node_flux_boundary(node, flux));
ek->ghost_communication();
{
auto const res = ek->get_node_is_boundary(node, true);
BOOST_REQUIRE(res);
Expand All @@ -186,6 +187,7 @@ BOOST_DATA_TEST_CASE(node_flux_boundary, bdata::make(all_eks()), ek_generator) {
}
{
BOOST_CHECK(ek->remove_node_from_flux_boundary(node));
ek->ghost_communication();
{
auto const res = ek->get_node_is_boundary(node, true);
BOOST_REQUIRE(res);
Expand All @@ -205,13 +207,16 @@ BOOST_DATA_TEST_CASE(node_flux_boundary, bdata::make(all_eks()), ek_generator) {
} else {
// Not in the local halo.
BOOST_CHECK(!ek->set_node_flux_boundary(node, flux));
ek->ghost_communication();
BOOST_CHECK(!ek->get_node_flux_at_boundary(node));
BOOST_CHECK(!ek->remove_node_from_flux_boundary(node));
ek->ghost_communication();
BOOST_CHECK(!ek->get_node_is_flux_boundary(node));
}
}

ek->clear_flux_boundaries();
ek->ghost_communication();
for (auto const &node : local_nodes_incl_ghosts(ek->get_lattice())) {
BOOST_CHECK(!(*ek->get_node_is_flux_boundary(node, true)));
}
Expand Down Expand Up @@ -321,6 +326,7 @@ BOOST_DATA_TEST_CASE(update_flux_boundary_from_shape, bdata::make(all_eks()),
std::vector<double> flux_flat(flux_3d.data(),
flux_3d.data() + flux_3d.num_elements());
ek->update_flux_boundary_from_shape(raster_flat, flux_flat);
ek->ghost_communication();
}

for (auto const &node : nodes) {
Expand Down
2 changes: 1 addition & 1 deletion testsuite/python/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -397,7 +397,7 @@ if(${ESPRESSO_TEST_NP} GREATER_EQUAL 6)
endif()
python_test(FILE ek_interface.py MAX_NUM_PROC 2)
python_test(FILE ek_diffusion.py MAX_NUM_PROC 1)
python_test(FILE ek_noflux.py MAX_NUM_PROC 1)
python_test(FILE ek_noflux.py MAX_NUM_PROC 2)
python_test(FILE ek_eof.py MAX_NUM_PROC 1)
python_test(FILE ek_fixedflux.py MAX_NUM_PROC 1)
python_test(FILE ek_bulk_reactions.py MAX_NUM_PROC 1)
Expand Down
4 changes: 2 additions & 2 deletions testsuite/python/ek_noflux.py
Original file line number Diff line number Diff line change
Expand Up @@ -28,7 +28,7 @@

@utx.skipIfMissingFeatures(["WALBERLA"])
class EKNoFlux(ut.TestCase):
BOX_L = 15.
BOX_L = 16.
AGRID = 1.0
DENSITY = 1
DIFFUSION_COEFFICIENT = 0.1
Expand Down Expand Up @@ -56,7 +56,7 @@ def detail_test_noflux(self, single_precision: bool):
decimal_precision: int = 7 if single_precision else 10

lattice = espressomd.electrokinetics.LatticeWalberla(
n_ghost_layers=1, agrid=self.AGRID)
n_ghost_layers=2, agrid=self.AGRID)

ekspecies = espressomd.electrokinetics.EKSpecies(
lattice=lattice, density=0.0, diffusion=self.DIFFUSION_COEFFICIENT,
Expand Down