25#include <blockforest/Initialization.h>
26#include <blockforest/StructuredBlockForest.h>
27#include <core/DataTypes.h>
28#include <core/cell/Cell.h>
29#include <domain_decomposition/IBlock.h>
34#include <initializer_list>
44 unsigned int n_ghost_layers)
45 : m_grid_dimensions{grid_dimensions}, m_node_grid{node_grid},
46 m_n_ghost_layers{n_ghost_layers} {
47 using walberla::real_t;
48 using walberla::uint_c;
50 for (
auto const i : {0u, 1u, 2u}) {
51 if (m_grid_dimensions[i] % node_grid[i] != 0) {
52 throw std::runtime_error(
53 "Lattice grid dimensions and MPI node grid are not compatible.");
55 if (m_grid_dimensions[i] % block_grid[i] != 0) {
56 throw std::runtime_error(
57 "Lattice grid dimensions and block grid are not compatible.");
61 auto constexpr lattice_constant = real_t{1};
62 auto const cells_per_block =
65 m_blocks = walberla::blockforest::createUniformBlockGrid(
67 uint_c(block_grid[0]), uint_c(block_grid[1]), uint_c(block_grid[2]),
69 uint_c(cells_per_block[0]), uint_c(cells_per_block[1]),
70 uint_c(cells_per_block[2]), lattice_constant,
72 uint_c(node_grid[0]), uint_c(node_grid[1]), uint_c(node_grid[2]),
77 for (IBlock &
block : *m_blocks) {
78 m_cached_blocks.push_back(&
block);
82[[nodiscard]] std::pair<Utils::Vector3d, Utils::Vector3d>
88 auto aa = to_vector3d(m_blocks->begin()->getAABB().min());
89 auto bb = to_vector3d(m_blocks->begin()->getAABB().max());
90 for (
auto const &
block : *m_blocks) {
91 auto cc =
block.getAABB();
92 for (
auto const i : {0u, 1u, 2u}) {
93 aa[i] = std::min(aa[i], cc.min()[i]);
94 bb[i] = std::max(bb[i], cc.max()[i]);
100[[nodiscard]] std::pair<Utils::Vector3i, Utils::Vector3i>
118 return ::walberla::get_block_and_cell(*
this, node,
false) != std::nullopt;
122 return ::walberla::get_block_and_cell(*
this, node,
true) != std::nullopt;
126 return ::walberla::get_block(*
this, pos,
false) !=
nullptr;
130 return ::walberla::get_block(*
this, pos,
true) !=
nullptr;
136 auto const grid_dimensions =
138 static_cast<int>(std::round(box_size[1] / agrid)),
139 static_cast<int>(std::round(box_size[2] / agrid))}};
140 for (
auto const i : {0u, 1u, 2u}) {
141 if (std::abs(grid_dimensions[i] * agrid - box_size[i]) / box_size[i] >
142 std::numeric_limits<double>::epsilon()) {
143 throw std::runtime_error(
144 "Box length not commensurate with agrid in direction " +
145 std::to_string(i) +
" length " + std::to_string(box_size[i]) +
146 " agrid " + std::to_string(agrid));
149 return grid_dimensions;
Vector implementation and trait types for boost qvm interoperability.
bool pos_in_local_halo(Utils::Vector3d const &pos) const
std::pair< Utils::Vector3i, Utils::Vector3i > get_local_grid_range() const
static Utils::Vector3i calc_grid_dimensions(Utils::Vector3d const &box_size, double agrid)
bool node_in_local_halo(Utils::Vector3i const &node) const
std::pair< Utils::Vector3d, Utils::Vector3d > get_local_domain() const
bool node_in_local_domain(Utils::Vector3i const &node) const
Utils::Vector3i get_block_corner(IBlock const &block, bool lower) const
LatticeWalberla(Utils::Vector3i const &grid_dimensions, Utils::Vector3i const &node_grid, Utils::Vector3i const &block_grid, unsigned int n_ghost_layers)
bool pos_in_local_domain(Utils::Vector3d const &pos) const
static double * block(double *p, std::size_t index, std::size_t size)
auto hadamard_division(Vector< T, N > const &a, Vector< U, N > const &b)
auto to_vector3d(Vector3< T > const &v) noexcept
auto convert_cell_corner_to_coord(real_vector auto const &corner)
Get the block-local coordinates of a block corner.
auto get_min_corner(IBlock const &block)
Get the block-local coordinates of the lower corner of a block.
auto get_max_corner(IBlock const &block)
Get the block-local coordinates of the upper corner of a block.