Loading [MathJax]/extensions/TeX/AMSmath.js
ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages Concepts
BlockAndCell.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2020-2025 The ESPResSo project
3 *
4 * This file is part of ESPResSo.
5 *
6 * ESPResSo is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation, either version 3 of the License, or
9 * (at your option) any later version.
10 *
11 * ESPResSo is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program. If not, see <http://www.gnu.org/licenses/>.
18 */
19
20#pragma once
21
22#include <blockforest/StructuredBlockForest.h>
23#include <core/DataTypes.h>
24#include <core/cell/Cell.h>
25#include <domain_decomposition/IBlock.h>
26
27#include "LatticeWalberla.hpp"
28
29#include <array>
30#include <cmath>
31#include <concepts>
32#include <memory>
33#include <optional>
34#include <type_traits>
35
36namespace detail {
37template <typename T> struct is_real_vector : std::false_type {};
38
39template <std::floating_point T>
40struct is_real_vector<std::array<T, 3>> : std::true_type {};
41
42template <std::floating_point T>
43struct is_real_vector<walberla::Vector3<T>> : std::true_type {};
44
45template <std::floating_point T>
46struct is_real_vector<Utils::Vector<T, 3>> : std::true_type {};
47} // namespace detail
48
49template <typename T>
50concept real_vector = detail::is_real_vector<T>::value;
51
52namespace walberla {
53
54inline Cell to_cell(Utils::Vector3i const &xyz) {
55 return {xyz[0], xyz[1], xyz[2]};
56}
57
59 IBlock *block;
61};
62
63template <typename T>
64IBlock *get_block_extended(LatticeWalberla const &lattice,
65 Utils::Vector<T, 3> const &pos,
66 unsigned int n_ghost_layers) {
67 auto const &cached_blocks = lattice.get_cached_blocks();
68 for (auto &block : cached_blocks) {
69 if (block->getAABB()
70 .getExtended(real_c(n_ghost_layers))
71 .contains(real_c(pos[0]), real_c(pos[1]), real_c(pos[2]))) {
72 return &(*block);
73 }
74 }
75 // Cell not in local blocks
76 return nullptr;
77}
78
79inline std::optional<BlockAndCell>
81 Utils::Vector3i const &node, bool consider_ghost_layers) {
82 auto const &blocks = lattice.get_blocks();
83 auto n_ghost_layers = 0u;
84 if (consider_ghost_layers) {
85 n_ghost_layers = lattice.get_ghost_layers();
86 }
87
88 auto block = get_block_extended(lattice, node, n_ghost_layers);
89 if (!block)
90 return std::nullopt;
91
92 // Transform coords to block local
93 Cell local_cell;
94
95 Cell global_cell = to_cell(node);
96 blocks->transformGlobalToBlockLocalCell(local_cell, *block, global_cell);
97 return {{block, local_cell}};
98}
99
100inline IBlock *get_block(::LatticeWalberla const &lattice,
101 Utils::Vector3d const &pos,
102 bool consider_ghost_layers) {
103 // Get block
104 auto const blocks = lattice.get_blocks();
105 auto block = blocks->getBlock(real_c(pos[0]), real_c(pos[1]), real_c(pos[2]));
106 if (consider_ghost_layers and !block) {
107 block = get_block_extended(lattice, pos, lattice.get_ghost_layers());
108 }
109 return block;
110}
111
112/**
113 * @brief Get the block-local coordinates of a block corner.
114 *
115 * This method leverages the fact that the grid spacing is unity in LB units,
116 * i.e. floating-point coordinates can be cast to integers indices.
117 */
118inline auto convert_cell_corner_to_coord(real_vector auto const &corner) {
119 return Utils::Vector3i{{static_cast<int>(std::round(corner[0])),
120 static_cast<int>(std::round(corner[1])),
121 static_cast<int>(std::round(corner[2]))}};
122}
123
124/** @brief Get the block-local coordinates of the lower corner of a block. */
125inline auto get_min_corner(IBlock const &block) {
126 return convert_cell_corner_to_coord(block.getAABB().minCorner());
127}
128
129/** @brief Get the block-local coordinates of the upper corner of a block. */
130inline auto get_max_corner(IBlock const &block) {
131 return convert_cell_corner_to_coord(block.getAABB().maxCorner());
132}
133
134} // namespace walberla
Definition Cell.hpp:96
Class that runs and controls the BlockForest in waLBerla.
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:172
Vector< T, 3 > Vector3
Definition Vector.hpp:161
\file PackInfoPdfDoublePrecision.cpp \author pystencils
std::optional< BlockAndCell > get_block_and_cell(::LatticeWalberla const &lattice, Utils::Vector3i const &node, bool consider_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_extended(LatticeWalberla const &lattice, Utils::Vector< T, 3 > const &pos, unsigned int n_ghost_layers)
auto get_min_corner(IBlock const &block)
Get the block-local coordinates of the lower corner of a block.
Cell to_cell(Utils::Vector3i const &xyz)
IBlock * get_block(::LatticeWalberla const &lattice, Utils::Vector3d const &pos, bool consider_ghost_layers)
auto get_max_corner(IBlock const &block)
Get the block-local coordinates of the upper corner of a block.