22#include <blockforest/StructuredBlockForest.h>
23#include <core/DataTypes.h>
24#include <core/cell/Cell.h>
25#include <domain_decomposition/IBlock.h>
29#include "../src/utils/types_conversion.hpp"
40template <
typename T>
struct is_real_vector : std::false_type {};
42template <std::
floating_po
int T>
43struct is_real_vector<
std::array<T, 3>> : std::true_type {};
45template <std::
floating_po
int T>
48template <std::
floating_po
int T>
49struct is_real_vector<
Utils::Vector<T, 3>> : std::true_type {};
56template <
typename T>
struct is_signed_integral_vector : std::false_type {};
58template <std::
signed_
integral T>
59struct is_signed_integral_vector<
std::array<T, 3>> : std::true_type {};
61template <std::
signed_
integral T>
62struct is_signed_integral_vector<
walberla::Vector3<T>> : std::true_type {};
64template <std::
signed_
integral T>
65struct is_signed_integral_vector<
Utils::Vector<T, 3>> : std::true_type {};
67template <>
struct is_signed_integral_vector<
walberla::
Cell> : std::true_type {
68 static_assert(std::integral<walberla::cell_idx_t> and
69 std::is_signed_v<walberla::cell_idx_t>);
79 return {xyz[0], xyz[1], xyz[2]};
88 unsigned int n_ghost_layers) {
90 for (
auto &
block : cached_blocks) {
92 .getExtended(real_c(n_ghost_layers))
93 .contains(real_c(pos[0]), real_c(pos[1]), real_c(pos[2]))) {
101inline std::optional<BlockAndCell>
104 bool consider_ghost_layers) {
106 auto n_ghost_layers = 0u;
107 if (consider_ghost_layers) {
119 blocks->transformGlobalToBlockLocalCell(local_cell, *
block, global_cell);
120 return {{
block, local_cell}};
125 bool consider_ghost_layers) {
128 auto block = blocks->getBlock(real_c(pos[0]), real_c(pos[1]), real_c(pos[2]));
129 if (consider_ghost_layers and !
block) {
143 static_cast<int>(std::round(corner[1])),
144 static_cast<int>(std::round(corner[2]))}};
157[[nodiscard]]
inline std::optional<walberla::cell::CellInterval>
161 auto const &cell_min = lower_corner;
165 if (not lower_bc or not upper_bc) {
169 auto const block_extent =
171 auto const global_lower_cell = lower_bc->cell;
172 auto const global_upper_cell = upper_bc->cell +
to_cell(block_extent);
173 return {CellInterval(global_lower_cell, global_upper_cell)};
177[[nodiscard]]
inline std::optional<walberla::cell::CellInterval>
183 if (not(upper_corner > block_lower_corner)) {
186 for (uint_t f = 0u; f < 3u; ++f) {
187 block_lower_corner[f] = std::max(block_lower_corner[f], lower_corner[f]);
190 if (not(block_upper_corner > lower_corner)) {
193 for (uint_t f = 0u; f < 3u; ++f) {
194 block_upper_corner[f] = std::min(block_upper_corner[f], upper_corner[f]);
197 auto const block_lower_cell =
to_cell(block_lower_corner - block_offset);
198 auto const block_upper_cell =
to_cell(block_upper_corner - block_offset);
199 return {CellInterval(block_lower_cell, block_upper_cell)};
216template <
typename Kernel>
220 auto const local_grid =
to_vector3i(ci.max() - ci.min() +
Cell(1, 1, 1));
221 auto const block_grid =
to_vector3i(bci.max() - bci.min() +
Cell(1, 1, 1));
222 auto const lower_cell = bci.min();
223 auto const upper_cell = bci.max();
228 for (
auto x = lower_cell.x(), i = 0; x <= upper_cell.x(); ++x, ++i) {
229 for (
auto y = lower_cell.y(), j = 0; y <= upper_cell.y(); ++y, ++j) {
230 for (
auto z = lower_cell.z(), k = 0; z <= upper_cell.z(); ++z, ++k) {
236 kernel(
static_cast<unsigned>(block_index),
237 static_cast<unsigned>(local_index), node);
Vector implementation and trait types for boost qvm interoperability.
Class that runs and controls the BlockForest in waLBerla.
Utils::Vector3i get_block_corner(IBlock const &block, bool lower) const
auto get_ghost_layers() const
auto const & get_cached_blocks() const
static DEVICE_QUALIFIER constexpr Vector< T, N > broadcast(typename Base::value_type const &value) noexcept
Create a vector that has all entries set to the same value.
static double * block(double *p, std::size_t index, std::size_t size)
int get_linear_index(int a, int b, int c, const Vector3i &adim)
\file PackInfoPdfDoublePrecision.cpp \author pystencils
IBlock * get_block_extended(LatticeWalberla const &lattice, auto const &pos, unsigned int n_ghost_layers)
auto convert_cell_corner_to_coord(real_vector auto const &corner)
Get the block-local coordinates of a block corner.
IBlock * get_block(::LatticeWalberla const &lattice, real_vector auto const &pos, bool consider_ghost_layers)
std::optional< BlockAndCell > get_block_and_cell(::LatticeWalberla const &lattice, signed_integral_vector auto const &node, bool consider_ghost_layers)
auto get_min_corner(IBlock const &block)
Get the block-local coordinates of the lower corner of a block.
Cell to_cell(signed_integral_vector auto const &xyz)
std::optional< walberla::cell::CellInterval > get_block_interval(::LatticeWalberla const &lattice, Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, Utils::Vector3i const &block_offset, IBlock const &block)
void copy_block_buffer(CellInterval const &bci, CellInterval const &ci, Utils::Vector3i const &block_offset, Utils::Vector3i const &lower_corner, Kernel &&kernel)
Synchronize data between a sliced block and a container.
Utils::Vector3i to_vector3i(Vector3< int > const &v) noexcept
std::optional< walberla::cell::CellInterval > get_interval(::LatticeWalberla const &lattice, Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner)
auto get_max_corner(IBlock const &block)
Get the block-local coordinates of the upper corner of a block.