ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
LBNodeAccess.impl.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2019-2026 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/**
23 * @file
24 * Out-of-class node access definitions for
25 * @ref walberla::LBWalberlaImpl.
26 */
27
28#include <utils/Vector.hpp>
29
30#include <array>
31#include <optional>
32#include <utility>
33#include <vector>
34
35namespace walberla {
36
37template <typename FloatType, lbmpy::Arch Architecture>
38std::optional<Utils::Vector3d>
40 Utils::Vector3i const &node, bool consider_ghosts) const {
41 assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::VEL)));
42 assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::UBB)));
43 if (m_has_boundaries) {
44 auto const is_boundary = get_node_is_boundary(node, consider_ghosts);
45 if (is_boundary and *is_boundary) {
46 return get_node_velocity_at_boundary(node, consider_ghosts);
47 }
48 }
49 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
50 if (!bc)
51 return std::nullopt;
52
53 auto field = bc->block->template uncheckedFastGetData<VectorField>(
54 m_velocity_field_id);
55 auto const vec = lbm::accessor::Vector::get(field, bc->cell);
56 return to_vector3d(vec);
57}
58
59template <typename FloatType, lbmpy::Arch Architecture>
61 Utils::Vector3i const &node, Utils::Vector3d const &v) {
62 m_pending_ghost_comm.set(GhostComm::PDF);
63 m_pending_ghost_comm.set(GhostComm::VEL);
64 auto bc = get_block_and_cell(get_lattice(), node, false);
65 if (!bc)
66 return false;
67
68 // We have to set both, the pdf and the stored velocity field
69 auto pdf_field = bc->block->template getData<PdfField>(m_pdf_field_id);
70 auto vel_field =
71 bc->block->template getData<VectorField>(m_velocity_field_id);
72 auto force_field =
73 bc->block->template getData<VectorField>(m_last_applied_force_field_id);
74 auto vel = to_vector3<FloatType>(v);
76 bc->cell);
77
78 return true;
79}
80
81template <typename FloatType, lbmpy::Arch Architecture>
83 Utils::Vector3i const &node, bool consider_ghosts) const {
84 assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::PDF)));
85 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
86 if (!bc)
87 return std::nullopt;
88
89 auto pdf_field =
90 bc->block->template uncheckedFastGetData<PdfField>(m_pdf_field_id);
91 auto const density =
92 lbm::accessor::Density::get(pdf_field, m_density, bc->cell);
93 return {double_c(density)};
94}
95
96template <typename FloatType, lbmpy::Arch Architecture>
98 Utils::Vector3i const &node, double density) {
99 m_pending_ghost_comm.set(GhostComm::PDF);
100 auto bc = get_block_and_cell(get_lattice(), node, false);
101 if (!bc)
102 return false;
103
104 auto pdf_field = bc->block->template getData<PdfField>(m_pdf_field_id);
105 lbm::accessor::Density::set(pdf_field, FloatType_c(density), m_density,
106 bc->cell);
107
108 return true;
109}
110
111template <typename FloatType, lbmpy::Arch Architecture>
112std::optional<std::vector<double>>
114 Utils::Vector3i const &node, bool consider_ghosts) const {
115 assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::PDF)));
116 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
117 if (!bc)
118 return std::nullopt;
119
120 auto pdf_field = bc->block->template getData<PdfField>(m_pdf_field_id);
121 auto const pop = lbm::accessor::Population::get(pdf_field, bc->cell);
122 std::vector<double> population(Stencil::Size);
123 for (uint_t f = 0u; f < Stencil::Size; ++f) {
124 population[f] = double_c(pop[f]);
125 }
126
127 return {std::move(population)};
128}
129
130template <typename FloatType, lbmpy::Arch Architecture>
132 Utils::Vector3i const &node, std::vector<double> const &population) {
133 m_pending_ghost_comm.set(GhostComm::PDF);
134 m_pending_ghost_comm.set(GhostComm::VEL);
135 auto bc = get_block_and_cell(get_lattice(), node, false);
136 if (!bc)
137 return false;
138
139 auto pdf_field = bc->block->template getData<PdfField>(m_pdf_field_id);
140 auto force_field =
141 bc->block->template getData<VectorField>(m_last_applied_force_field_id);
142 auto vel_field =
143 bc->block->template getData<VectorField>(m_velocity_field_id);
144 std::array<FloatType, Stencil::Size> pop;
145 for (uint_t f = 0u; f < Stencil::Size; ++f) {
146 pop[f] = FloatType_c(population[f]);
147 }
149 bc->cell);
150
151 return true;
152}
153
154template <typename FloatType, lbmpy::Arch Architecture>
155std::optional<Utils::Vector3d>
157 Utils::Vector3i const &node) const {
158 auto const bc = get_block_and_cell(get_lattice(), node, true);
159 if (!bc)
160 return std::nullopt;
161
162 auto field =
163 bc->block->template getData<VectorField>(m_force_to_be_applied_id);
164 auto const vec = lbm::accessor::Vector::get(field, bc->cell);
165 return zero_centered_to_md(to_vector3d(vec));
166}
167
168template <typename FloatType, lbmpy::Arch Architecture>
169std::optional<Utils::Vector3d>
171 Utils::Vector3i const &node, bool consider_ghosts) const {
172 assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::LAF)));
173 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
174 if (!bc)
175 return std::nullopt;
176
177 auto const field =
178 bc->block->template getData<VectorField>(m_last_applied_force_field_id);
179 auto const vec = lbm::accessor::Vector::get(field, bc->cell);
180 return zero_centered_to_md(to_vector3d(vec));
181}
182
183template <typename FloatType, lbmpy::Arch Architecture>
185 Utils::Vector3i const &node, Utils::Vector3d const &force) {
186 m_pending_ghost_comm.set(GhostComm::VEL);
187 m_pending_ghost_comm.set(GhostComm::LAF);
188 auto bc = get_block_and_cell(get_lattice(), node, false);
189 if (!bc)
190 return false;
191
192 auto pdf_field = bc->block->template getData<PdfField>(m_pdf_field_id);
193 auto force_field =
194 bc->block->template getData<VectorField>(m_last_applied_force_field_id);
195 auto vel_field =
196 bc->block->template getData<VectorField>(m_velocity_field_id);
197 auto const vec = to_vector3<FloatType>(force);
199 bc->cell);
200
201 return true;
202}
203
204template <typename FloatType, lbmpy::Arch Architecture>
205std::optional<Utils::VectorXd<9>>
207 Utils::Vector3i const &node) const {
208 auto bc = get_block_and_cell(get_lattice(), node, false);
209 if (!bc)
210 return std::nullopt;
211
212 auto pdf_field = bc->block->template getData<PdfField>(m_pdf_field_id);
213 auto tensor =
215 pressure_tensor_correction(tensor);
216 return to_vector9d(tensor);
217}
218
219} // namespace walberla
Vector implementation and trait types for boost qvm interoperability.
std::optional< Utils::Vector3d > get_node_last_applied_force(Utils::Vector3i const &node, bool consider_ghosts=false) const override
std::optional< Utils::Vector3d > get_node_velocity(Utils::Vector3i const &node, bool consider_ghosts=false) const override
bool set_node_last_applied_force(Utils::Vector3i const &node, Utils::Vector3d const &force) override
bool set_node_density(Utils::Vector3i const &node, double density) override
std::optional< std::vector< double > > get_node_population(Utils::Vector3i const &node, bool consider_ghosts=false) const override
std::optional< double > get_node_density(Utils::Vector3i const &node, bool consider_ghosts=false) const override
std::optional< Utils::Vector3d > get_node_force_to_be_applied(Utils::Vector3i const &node) const override
bool set_node_velocity(Utils::Vector3i const &node, Utils::Vector3d const &v) override
bool set_node_population(Utils::Vector3i const &node, std::vector< double > const &population) override
std::optional< Utils::VectorXd< 9 > > get_node_pressure_tensor(Utils::Vector3i const &node) const override
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.
void set(GhostLayerField< double, uint_t{19u}> *pdf_field, double const rho_in, double const density, Cell const &cell)
double get(GhostLayerField< double, uint_t{19u}> const *pdf_field, double const density, Cell const &cell)
void set(GhostLayerField< double, uint_t{19u}> const *pdf_field, GhostLayerField< double, uint_t{3u}> *velocity_field, GhostLayerField< double, uint_t{3u}> *force_field, Vector3< double > const &force, double const density, Cell const &cell)
void set(GhostLayerField< double, uint_t{19u}> *pdf_field, std::array< double, 19u > const &pop, Cell const &cell)
auto get(GhostLayerField< double, uint_t{19u}> const *pdf_field, Cell const &cell)
auto get(GhostLayerField< double, uint_t{19u}> const *pdf_field, double const density, Cell const &cell)
auto get(GhostLayerField< double, uint_t{3u}> const *vec_field, Cell const &cell)
void set(GhostLayerField< double, uint_t{19u}> *pdf_field, GhostLayerField< double, uint_t{3u}> *velocity_field, GhostLayerField< double, uint_t{3u}> const *force_field, Vector3< double > const &u, Cell const &cell)
\file PackInfoPdfDoublePrecision.cpp \author pystencils
auto to_vector3d(Vector3< T > const &v) noexcept
std::optional< BlockAndCell > get_block_and_cell(::LatticeWalberla const &lattice, signed_integral_vector auto const &node, bool consider_ghost_layers)
auto to_vector9d(Matrix3< T > const &m) noexcept