ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
BlockAndCell.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2020-2023 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 <memory>
30#include <optional>
31
32namespace walberla {
33// Helpers to retrieve blocks and cells
35 IBlock *block;
37};
38
39template <typename T>
40IBlock *get_block_extended(LatticeWalberla const &lattice,
41 Utils::Vector<T, 3> const &pos,
42 unsigned int n_ghost_layers) {
43 auto const &cached_blocks = lattice.get_cached_blocks();
44 for (auto &block : cached_blocks) {
45 if (block->getAABB()
46 .getExtended(real_c(n_ghost_layers))
47 .contains(real_c(pos[0]), real_c(pos[1]), real_c(pos[2]))) {
48 return &(*block);
49 }
50 }
51 // Cell not in local blocks
52 return nullptr;
53}
54
55inline std::optional<BlockAndCell>
57 Utils::Vector3i const &node, bool consider_ghost_layers) {
58 auto const &blocks = lattice.get_blocks();
59 auto n_ghost_layers = 0u;
60 if (consider_ghost_layers) {
61 n_ghost_layers = lattice.get_ghost_layers();
62 }
63
64 auto block = get_block_extended(lattice, node, n_ghost_layers);
65 if (!block)
66 return std::nullopt;
67
68 // Transform coords to block local
69 Cell local_cell;
70
71 Cell global_cell{uint_c(node[0]), uint_c(node[1]), uint_c(node[2])};
72 blocks->transformGlobalToBlockLocalCell(local_cell, *block, global_cell);
73 return {{block, local_cell}};
74}
75
76inline IBlock *get_block(::LatticeWalberla const &lattice,
77 Utils::Vector3d const &pos,
78 bool consider_ghost_layers) {
79 // Get block
80 auto const blocks = lattice.get_blocks();
81 auto block = blocks->getBlock(real_c(pos[0]), real_c(pos[1]), real_c(pos[2]));
82 if (consider_ghost_layers and !block) {
83 block = get_block_extended(lattice, pos, lattice.get_ghost_layers());
84 }
85 return block;
86}
87
88} // namespace walberla
Definition Cell.hpp:97
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
\file PackInfoPdfDoublePrecision.cpp \author pystencils
std::optional< BlockAndCell > get_block_and_cell(::LatticeWalberla const &lattice, Utils::Vector3i const &node, bool consider_ghost_layers)
IBlock * get_block_extended(LatticeWalberla const &lattice, Utils::Vector< T, 3 > const &pos, unsigned int n_ghost_layers)
IBlock * get_block(::LatticeWalberla const &lattice, Utils::Vector3d const &pos, bool consider_ghost_layers)