ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
LBWalberlaImpl.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 * @ref walberla::LBWalberlaImpl implements the interface of the LB
25 * waLBerla bridge using sweeps generated by lbmpy
26 * (see <tt>maintainer/walberla_kernels</tt>).
27 */
28
29#include <blockforest/Initialization.h>
30#include <blockforest/StructuredBlockForest.h>
31#include <domain_decomposition/BlockDataID.h>
32#include <domain_decomposition/IBlock.h>
33#include <field/AddToStorage.h>
34#include <field/vtk/FlagFieldCellFilter.h>
35#include <field/vtk/VTKWriter.h>
36#include <stencil/D3Q19.h>
37#include <stencil/D3Q27.h>
38#include <waLBerlaDefinitions.h>
39#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
40#include <gpu/AddGPUFieldToStorage.h>
41#endif
42
43#include "../BoundaryHandling.hpp"
44#include "../BoundaryPackInfo.hpp"
45#include "../utils/boundary.hpp"
46#include "../utils/types_conversion.hpp"
48#include "ResetForce.hpp"
49#include "lb_fields.hpp"
50#include "lb_kernels.hpp"
51#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
52#include "lb_fields.cuh"
53#include "lb_kernels.cuh"
54#endif
55
63
64#include <utils/Vector.hpp>
65
66#include <array>
67#include <bitset>
68#include <cstddef>
69#include <functional>
70#include <initializer_list>
71#include <limits>
72#include <memory>
73#include <optional>
74#include <stdexcept>
75#include <string>
76#include <type_traits>
77#include <utility>
78#include <variant>
79#include <vector>
80
81namespace walberla {
82
83/** @brief Class that runs and controls the LB on waLBerla. */
84template <typename FloatType, lbmpy::Arch Architecture>
86#if not defined(WALBERLA_BUILD_WITH_CUDA)
87 static_assert(Architecture != lbmpy::Arch::GPU,
88 "waLBerla was compiled without CUDA support");
89#endif
90protected:
91 // ---- Types & Constants ----
92
93 using Kernels = detail::KernelTrait<FloatType, Architecture>;
95 typename Kernels::DynamicUBB>;
97 std::variant<typename Kernels::StreamCollisionModelThermalized,
99
100public:
101 /** @brief Stencil for collision and streaming operations. */
102 using Stencil = stencil::D3Q19;
103 /** @brief Stencil for ghost communication (includes domain corners). */
104 using StencilFull = stencil::D3Q27;
105 /** @brief Lattice model (e.g. blockforest). */
107
108protected:
109 // "underlying" field types (`GPUField` has no f-size info at compile time)
112
113public:
117#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
118 using GPUField = gpu::GPUField<FloatType>;
119#endif
120
121 struct GhostComm {
122 /** @brief Ghost communication operations. */
123 enum GhostCommFlags : unsigned {
124 PDF, ///< PDFs communication
125 VEL, ///< velocities communication
126 LAF, ///< last applied forces communication
127 UBB, ///< boundaries communication
128 SIZE
129 };
130 };
131
132protected:
133 /**
134 * @brief Full communicator.
135 * We use the D3Q27 directions to update cells along the diagonals during
136 * a full ghost communication. This is needed to properly update the corners
137 * of the ghost layer when setting cell velocities or populations.
138 */
140 FieldTrait<FloatType, Stencil,
141 Architecture>::template RegularCommScheme<stencil::D3Q27>;
143 FieldTrait<FloatType, Stencil,
144 Architecture>::template BoundaryCommScheme<stencil::D3Q27>;
145 /**
146 * @brief Regular communicator.
147 * We use the same directions as the stencil during integration.
148 */
150 FieldTrait<FloatType, Stencil,
151 Architecture>::template RegularCommScheme<Stencil>;
152 template <class Field>
153 using PackInfo =
155
156protected:
157 // ---- Member Variables ----
158
159 // Physical parameters
160 FloatType m_viscosity; /// kinematic viscosity
161 FloatType m_density;
162 FloatType m_kT;
163 unsigned int m_seed;
164 double m_zc_to_md; // zero-centered conversion factor to MD units
165 double m_zc_to_lb; // zero-centered conversion factor to LB units
166
167 // lattice
168 std::shared_ptr<LatticeWalberla> m_lattice;
169
170 // Block data access handles
171 BlockDataID m_pdf_field_id;
173 BlockDataID m_flag_field_id;
174
177
180
181 std::optional<BlockDataID> m_pressure_tensor_field_id;
182#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
183 std::optional<BlockDataID> m_pdf_cpu_field_id;
184 std::optional<BlockDataID> m_vel_cpu_field_id;
185#endif
186
187 /** Flag for boundary cells. */
188 FlagUID const Boundary_flag{"boundary"};
189 bool m_has_boundaries{false};
190
191 // boundaries
192 std::shared_ptr<BoundaryModel> m_boundary;
193
194 // communicators
195 std::shared_ptr<BoundaryFullCommunicator> m_boundary_communicator;
196 std::shared_ptr<RegularFullCommunicator> m_full_communicator;
197 std::shared_ptr<RegularFullCommunicator> m_pdf_communicator;
198 std::shared_ptr<RegularFullCommunicator> m_vel_communicator;
199 std::shared_ptr<RegularFullCommunicator> m_laf_communicator;
200 std::shared_ptr<PDFStreamingCommunicator> m_pdf_streaming_communicator;
201 std::bitset<GhostComm::SIZE> m_pending_ghost_comm;
203
204 // collision sweep
205 std::shared_ptr<CollisionModel> m_collision_model;
206
207 // force reset sweep + external force handling
208 std::shared_ptr<ResetForce<PdfField, VectorField>> m_reset_force;
209
210 // velocity update sweep
211 std::shared_ptr<typename Kernels::UpdateVelFromPDF>
213
214 // Lees-Edwards boundary interpolation
215 std::shared_ptr<LeesEdwardsPack> m_lees_edwards_callbacks;
216 std::shared_ptr<InterpolateAndShiftAtBoundary<_PdfField, FloatType>>
218 std::shared_ptr<InterpolateAndShiftAtBoundary<_VectorField, FloatType>>
220 std::shared_ptr<InterpolateAndShiftAtBoundary<_VectorField, FloatType>>
222
223public:
224 template <typename T> FloatType FloatType_c(T t) const {
225 return numeric_cast<FloatType>(t);
226 }
227
228 [[nodiscard]] std::size_t stencil_size() const noexcept override {
229 return static_cast<std::size_t>(Stencil::Size);
230 }
231
232 [[nodiscard]] bool is_double_precision() const noexcept override {
233 return std::is_same_v<FloatType, double>;
234 }
235
236 [[nodiscard]] bool is_gpu() const noexcept override {
237 return Architecture == lbmpy::Arch::GPU;
238 }
239
240public:
241 LBWalberlaImpl(std::shared_ptr<LatticeWalberla> lattice, double viscosity,
242 double density)
244 m_kT(FloatType{0}), m_seed(0u), m_zc_to_md(density),
245 m_zc_to_lb(1. / density), m_lattice(std::move(lattice)),
247
248 auto const &blocks = m_lattice->get_blocks();
249 auto const n_ghost_layers = m_lattice->get_ghost_layers();
250 if (n_ghost_layers == 0u)
251 throw std::runtime_error("At least one ghost layer must be used");
252
253 // Initialize and register fields (must use the "underlying" types)
254 m_pdf_field_id = add_to_storage<_PdfField>("pdfs");
255 m_pdf_tmp_field_id = add_to_storage<_PdfField>("pdfs_tmp");
256 m_last_applied_force_field_id = add_to_storage<_VectorField>("force last");
257 m_force_to_be_applied_id = add_to_storage<_VectorField>("force next");
258 m_velocity_field_id = add_to_storage<_VectorField>("velocity");
259 m_vel_tmp_field_id = add_to_storage<_VectorField>("velocity_tmp");
260
261 // Initialize and register pdf field with zero centered density
262 auto pdf_setter = typename Kernels::InitialPDFsSetter(
264 for (auto &block : *blocks) {
265 pdf_setter(&block);
266 }
267
268 // Initialize and register flag field (fluid/boundary)
269 m_flag_field_id = field::addFlagFieldToStorage<FlagField>(
270 blocks, "flag field", n_ghost_layers);
271 // Initialize boundary sweep
272 reset_boundary_handling(m_lattice->get_blocks());
273
274 // Set up the communication and register fields
277
279
280 // Instantiate the sweep responsible for force double buffering and
281 // external forces
282 m_reset_force = std::make_shared<ResetForce<PdfField, VectorField>>(
284
285 // Instantiate velocity update sweep
287 std::make_shared<typename Kernels::UpdateVelFromPDF>(
289 }
290
291 ~LBWalberlaImpl() override = default;
292
293 // ---- Integration (Core LB Algorithm) ----
294
295 void integrate() override {
296 integrate_pull_scheme();
298 }
299
300protected:
301 void integrate_vtk_writers() override {
302 for (auto const &it : m_vtk_auto) {
303 auto &vtk_handle = it.second;
304 if (vtk_handle->enabled) {
305 vtk::writeFiles(vtk_handle->ptr)();
306 vtk_handle->execution_count++;
307 }
308 }
309 }
310
311private:
312 /**
313 * @brief One LB time step using the pull scheme.
314 * Sequence: reset forces, stream-collide, communicate PDFs,
315 * apply Lees-Edwards interpolation (if active), handle boundaries,
316 * update velocity field from PDFs.
317 */
318 void integrate_pull_scheme() {
320 auto const &blocks = get_lattice().get_blocks();
321 // Reset force fields
322 integrate_reset_force(blocks);
323 // LB stream collide
324 integrate_stream_collide(blocks);
325 // Mark pending ghost layer updates
326 // As pdf and laf are communicated directly afterwards, they are not set
329 m_pdf_streaming_communicator->communicate();
330 if (has_lees_edwards_bc()) {
331 apply_lees_edwards_pdf_interpolation(blocks);
332 apply_lees_edwards_last_applied_force_interpolation(blocks);
333 }
334 // Handle boundaries
335 if (m_has_boundaries) {
336 integrate_boundaries(blocks);
337 }
338 // Update velocities from pdfs
339 integrate_update_velocities_from_pdf(blocks);
340
341 if (has_lees_edwards_bc()) {
342 apply_lees_edwards_vel_interpolation_and_shift(blocks);
343 }
344 }
345
346 void integrate_stream_collide(std::shared_ptr<BlockStorage> const &blocks) {
347 auto &cm_variant = *m_collision_model;
348 for (auto &block : *blocks) {
349 auto const block_variant = std::variant<IBlock *>(&block);
350 std::visit(m_run_stream_collide_sweep, cm_variant, block_variant);
351 }
352 if (auto *cm =
353 std::get_if<typename Kernels::StreamCollisionModelThermalized>(
354 &cm_variant)) {
355 cm->setTime_step(cm->getTime_step() + 1u);
356 }
357 }
358
359 void integrate_reset_force(std::shared_ptr<BlockStorage> const &blocks) {
360 for (auto &block : *blocks)
361 (*m_reset_force)(&block);
362 }
363
364 void integrate_boundaries(std::shared_ptr<BlockStorage> const &blocks) {
365 for (auto &block : *blocks)
366 (*m_boundary)(&block);
367 }
368
369 void integrate_update_velocities_from_pdf(
370 std::shared_ptr<BlockStorage> const &blocks) {
371 for (auto &block : *blocks)
373 }
374
375private:
376 // ---- Collision Model ----
377
378 /**
379 * @brief Visitor for dispatching stream-collide sweeps.
380 * Handles both thermalized and Lees-Edwards collision models
381 * via @c std::visit on the @ref CollisionModel variant.
382 */
383 class StreamCollideSweepVisitor {
384 public:
385 using StructuredBlockStorage = LatticeWalberla::Lattice_T;
386
387 void operator()(typename Kernels::StreamCollisionModelThermalized &cm,
388 IBlock *b) {
389 cm.configure(m_storage, b);
390 cm(b);
391 }
392
393 void operator()(typename Kernels::StreamCollisionModelLeesEdwards &cm,
394 IBlock *b) {
395 cm.setV_s(static_cast<decltype(cm.getV_s())>(
396 m_lees_edwards_callbacks->get_shear_velocity()));
397 cm(b);
398 }
399
400 StreamCollideSweepVisitor() = default;
401 StreamCollideSweepVisitor(std::shared_ptr<StructuredBlockStorage> storage) {
402 m_storage = std::move(storage);
403 }
404 StreamCollideSweepVisitor(std::shared_ptr<StructuredBlockStorage> storage,
405 std::shared_ptr<LeesEdwardsPack> callbacks) {
406 m_storage = std::move(storage);
407 m_lees_edwards_callbacks = std::move(callbacks);
408 }
409
410 private:
411 std::shared_ptr<StructuredBlockStorage> m_storage{};
412 std::shared_ptr<LeesEdwardsPack> m_lees_edwards_callbacks{};
413 };
414 StreamCollideSweepVisitor m_run_stream_collide_sweep{};
415
416 /** @brief Relaxation rate omega from kinematic viscosity: 2/(6*nu+1). */
417 FloatType shear_mode_relaxation_rate() const;
418 /**
419 * @brief Odd-mode relaxation rate for the magic parameter relation.
420 * Ensures optimal bounce-back wall location for the two-relaxation-time
421 * model. Default magic number is 3/16.
422 */
423 FloatType odd_mode_relaxation_rate(
424 FloatType shear_relaxation,
425 FloatType magic_number = FloatType{3} / FloatType{16}) const;
426
427public:
428 void set_collision_model(double kT, unsigned int seed) override;
430 std::unique_ptr<LeesEdwardsPack> &&lees_edwards_pack) override;
431 void check_lebc(unsigned int shear_direction,
432 unsigned int shear_plane_normal) const override;
433
434public:
435 // ---- Ghost Communication ----
436
437 /**
438 * @brief Perform all pending ghost layer updates.
439 * Uses a lazy scheme: ghost communications are only executed when
440 * they have been marked as pending by a preceding write operation.
441 */
451
452 void ghost_communication_pdf() override {
455 m_pdf_communicator->communicate();
456 if (has_lees_edwards_bc()) {
457 auto const &blocks = get_lattice().get_blocks();
458 apply_lees_edwards_pdf_interpolation(blocks);
459 }
461 }
462 }
463
464 void ghost_communication_vel() override {
467 m_vel_communicator->communicate();
468 if (has_lees_edwards_bc()) {
469 auto const &blocks = get_lattice().get_blocks();
470 apply_lees_edwards_vel_interpolation_and_shift(blocks);
471 }
473 }
474 }
475
476 void ghost_communication_laf() override {
479 m_laf_communicator->communicate();
480 if (has_lees_edwards_bc()) {
481 auto const &blocks = get_lattice().get_blocks();
482 apply_lees_edwards_last_applied_force_interpolation(blocks);
483 }
485 }
486 }
487
495
496 /** @brief Communicate all fields at once using the D3Q27 stencil. */
499 m_full_communicator->communicate();
500 if (has_lees_edwards_bc()) {
502 }
506 }
507
508private:
509 // ---- Lees-Edwards Boundary Conditions ----
510
511 auto has_lees_edwards_bc() const {
512 return std::holds_alternative<
513 typename Kernels::StreamCollisionModelLeesEdwards>(*m_collision_model);
514 }
515
516 void apply_lees_edwards_pdf_interpolation(
517 std::shared_ptr<BlockStorage> const &blocks) {
518 for (auto &block : *blocks)
520 }
521
522 void apply_lees_edwards_vel_interpolation_and_shift(
523 std::shared_ptr<BlockStorage> const &blocks) {
524 for (auto &block : *blocks)
526 }
527
528 void apply_lees_edwards_last_applied_force_interpolation(
529 std::shared_ptr<BlockStorage> const &blocks) {
530 for (auto &block : *blocks)
532 }
533
534public:
536 auto const &blocks = get_lattice().get_blocks();
537 apply_lees_edwards_pdf_interpolation(blocks);
538 apply_lees_edwards_vel_interpolation_and_shift(blocks);
539 apply_lees_edwards_last_applied_force_interpolation(blocks);
540 }
541
542public:
543 // ---- Node & Slice Accessors (by quantity) ----
544
545 // Velocity
546 std::optional<Utils::Vector3d>
548 bool consider_ghosts = false) const override;
549 bool set_node_velocity(Utils::Vector3i const &node,
550 Utils::Vector3d const &v) override;
551 std::vector<double>
552 get_slice_velocity(Utils::Vector3i const &lower_corner,
553 Utils::Vector3i const &upper_corner) const override;
554 void set_slice_velocity(Utils::Vector3i const &lower_corner,
555 Utils::Vector3i const &upper_corner,
556 std::vector<double> const &velocity) override;
557
558 // Density
559 std::optional<double>
561 bool consider_ghosts = false) const override;
562 bool set_node_density(Utils::Vector3i const &node, double density) override;
563 std::vector<double>
564 get_slice_density(Utils::Vector3i const &lower_corner,
565 Utils::Vector3i const &upper_corner) const override;
566 void set_slice_density(Utils::Vector3i const &lower_corner,
567 Utils::Vector3i const &upper_corner,
568 std::vector<double> const &density) override;
569
570 // Population
571 std::optional<std::vector<double>>
573 bool consider_ghosts = false) const override;
574 bool set_node_population(Utils::Vector3i const &node,
575 std::vector<double> const &population) override;
576 std::vector<double>
577 get_slice_population(Utils::Vector3i const &lower_corner,
578 Utils::Vector3i const &upper_corner) const override;
579 void set_slice_population(Utils::Vector3i const &lower_corner,
580 Utils::Vector3i const &upper_corner,
581 std::vector<double> const &population) override;
582
583 // Force
584 std::optional<Utils::Vector3d>
585 get_node_force_to_be_applied(Utils::Vector3i const &node) const override;
586 std::optional<Utils::Vector3d>
588 bool consider_ghosts = false) const override;
590 Utils::Vector3d const &force) override;
591 std::vector<double> get_slice_last_applied_force(
592 Utils::Vector3i const &lower_corner,
593 Utils::Vector3i const &upper_corner) const override;
594 void set_slice_last_applied_force(Utils::Vector3i const &lower_corner,
595 Utils::Vector3i const &upper_corner,
596 std::vector<double> const &force) override;
597
598 // Pressure tensor
599 std::optional<Utils::VectorXd<9>>
600 get_node_pressure_tensor(Utils::Vector3i const &node) const override;
601 std::vector<double>
602 get_slice_pressure_tensor(Utils::Vector3i const &lower_corner,
603 Utils::Vector3i const &upper_corner) const override;
604
605private:
606 // ---- Interpolation (position-based access) ----
607
608 /** @brief Return a B-spline interpolation kernel for force distribution. */
609 auto make_force_interpolation_kernel() const;
610 /** @brief Return a B-spline interpolation kernel for velocity readout. */
611 auto make_velocity_interpolation_kernel() const;
612 /** @brief Return a B-spline interpolation kernel for density readout. */
613 auto make_density_interpolation_kernel() const;
614
615public:
616 std::function<bool(Utils::Vector3d const &)>
617 make_lattice_position_checker(bool consider_points_in_halo) const override;
618 bool add_force_at_pos(Utils::Vector3d const &pos,
619 Utils::Vector3d const &force) override;
620 void add_forces_at_pos(std::vector<Utils::Vector3d> const &pos,
621 std::vector<Utils::Vector3d> const &forces) override;
622 std::optional<Utils::Vector3d>
624 bool consider_points_in_halo = false) const override;
625 std::vector<Utils::Vector3d>
626 get_velocities_at_pos(std::vector<Utils::Vector3d> const &pos) override;
627 std::optional<double>
629 bool consider_points_in_halo = false) const override;
630 std::vector<double>
631 get_densities_at_pos(std::vector<Utils::Vector3d> const &pos) override;
632
633public:
634 // ---- Boundary Handling ----
635
636 void reset_boundary_handling(std::shared_ptr<BlockStorage> const &blocks) {
637 auto const [lc, uc] = m_lattice->get_local_grid_range(true);
638 m_boundary =
639 std::make_shared<BoundaryModel>(blocks, m_pdf_field_id, m_flag_field_id,
640 CellInterval{to_cell(lc), to_cell(uc)});
641 }
642
643 void on_boundary_add();
644 void clear_boundaries() override;
645 void reallocate_ubb_field() override;
646 void
647 update_boundary_from_shape(std::vector<int> const &raster_flat,
648 std::vector<double> const &data_flat) override;
649 std::optional<Utils::Vector3d>
651 bool consider_ghosts = false) const override;
653 Utils::Vector3d const &velocity) override;
654 std::vector<std::optional<Utils::Vector3d>> get_slice_velocity_at_boundary(
655 Utils::Vector3i const &lower_corner,
656 Utils::Vector3i const &upper_corner) const override;
658 Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
659 std::vector<std::optional<Utils::Vector3d>> const &velocity) override;
660 std::optional<Utils::Vector3d>
661 get_node_boundary_force(Utils::Vector3i const &node) const override;
662 bool remove_node_from_boundary(Utils::Vector3i const &node) override;
663 std::optional<bool>
665 bool consider_ghosts = false) const override;
666 std::vector<bool>
667 get_slice_is_boundary(Utils::Vector3i const &lower_corner,
668 Utils::Vector3i const &upper_corner) const override;
670 std::vector<int> const &raster_flat) const override;
671 [[nodiscard]] Utils::Vector3d get_boundary_force() const override;
672
673private:
674 [[nodiscard]] Utils::Vector3i flat_index_to_node(int index) const;
675 [[nodiscard]] Utils::Vector3i get_neighbor_node(Utils::Vector3i const &node,
676 int dir) const;
677
678public:
679 // ---- Global Reductions & Physical Parameters ----
680
681 // Global pressure tensor
682 [[nodiscard]] Utils::VectorXd<9> get_pressure_tensor() const override {
683 Matrix3<FloatType> tensor(FloatType{0});
684 for (auto const &block : *get_lattice().get_blocks()) {
685 auto pdf_field = block.template getData<PdfField>(m_pdf_field_id);
687 }
688 auto const &grid_size = get_lattice().get_grid_dimensions();
689 auto const number_of_nodes = Utils::product(grid_size);
691 return to_vector9d(tensor) * (1. / static_cast<double>(number_of_nodes));
692 }
693
694 // Global momentum
695 [[nodiscard]] Utils::Vector3d get_momentum() const override {
696 Vector3<FloatType> mom(FloatType{0});
697 for (auto const &block : *get_lattice().get_blocks()) {
698 auto pdf_field = block.template getData<PdfField>(m_pdf_field_id);
699 auto force_field =
700 block.template getData<VectorField>(m_last_applied_force_field_id);
701 mom += lbm::accessor::MomentumDensity::reduce(pdf_field, force_field,
702 m_density);
703 }
704 return to_vector3d(mom);
705 }
706
707 // Global external force
708 void set_external_force(Utils::Vector3d const &ext_force) override {
709 m_reset_force->set_ext_force(zero_centered_to_lb(ext_force));
710 }
711
712 [[nodiscard]] Utils::Vector3d get_external_force() const noexcept override {
713 return zero_centered_to_md(m_reset_force->get_ext_force());
714 }
715
716 void set_viscosity(double viscosity) override {
717 m_viscosity = FloatType_c(viscosity);
718 }
719
720 [[nodiscard]] double get_viscosity() const noexcept override {
721 return static_cast<double>(m_viscosity);
722 }
723
724 [[nodiscard]] double get_density() const noexcept override {
725 return static_cast<double>(m_density);
726 }
727
728 [[nodiscard]] double get_kT() const noexcept override {
729 return static_cast<double>(m_kT);
730 }
731
732 [[nodiscard]] unsigned int get_seed() const noexcept override {
733 return m_seed;
734 }
735
736 [[nodiscard]] std::optional<uint64_t> get_rng_state() const override {
737 auto const cm =
738 std::get_if<typename Kernels::StreamCollisionModelThermalized>(
740 if (!cm or m_kT == 0.) {
741 return std::nullopt;
742 }
743 return {static_cast<uint64_t>(cm->getTime_step())};
744 }
745
746 void set_rng_state(uint64_t counter) override {
747 auto const cm =
748 std::get_if<typename Kernels::StreamCollisionModelThermalized>(
750 if (!cm or m_kT == 0.) {
751 throw std::runtime_error("This LB instance is unthermalized");
752 }
753 assert(counter <=
754 static_cast<uint32_t>(std::numeric_limits<uint_t>::max()));
755 cm->setTime_step(static_cast<uint32_t>(counter));
756 }
757
758 [[nodiscard]] LatticeWalberla const &get_lattice() const noexcept override {
759 return *m_lattice;
760 }
761
762 [[nodiscard]] std::size_t get_velocity_field_id() const noexcept override {
763 return m_velocity_field_id;
764 }
765
766 [[nodiscard]] std::size_t get_force_field_id() const noexcept override {
768 }
769
770 /**
771 * @brief Correction factor for off-diagonal pressure tensor elements.
772 * Compensates for the viscosity-dependent error in the non-equilibrium
773 * stress: factor = nu / (nu + 1/6).
774 */
776 return m_viscosity / (m_viscosity + FloatType{1} / FloatType{6});
777 }
778
779 void pressure_tensor_correction(Matrix3<FloatType> &tensor) const {
780 auto const revert_factor = pressure_tensor_correction_factor();
781 for (auto const i : {1u, 2u, 3u, 5u, 6u, 7u}) {
782 tensor[i] *= revert_factor;
783 }
784 }
785
786 void pressure_tensor_correction(std::span<FloatType, 9ul> tensor) const {
787 auto const revert_factor = pressure_tensor_correction_factor();
788 for (auto const i : {1u, 2u, 3u, 5u, 6u, 7u}) {
789 tensor[i] *= revert_factor;
790 }
791 }
792
793protected:
794 /**
795 * @brief Scale data by a conversion factor (in-place).
796 * Used for zero-centered density representation: LB internally stores
797 * density fluctuations around zero, while the user interface uses
798 * absolute densities. The conversion factors @ref m_zc_to_md and
799 * @ref m_zc_to_lb translate between these representations.
800 */
801 template <typename T>
802 void zero_centered_transform_impl(T &data, auto const factor) const {
803 if constexpr (std::is_arithmetic_v<T>) {
804 static_assert(std::is_floating_point_v<T>);
805 data *= static_cast<T>(factor);
806 } else {
807 auto const coef = static_cast<typename T::value_type>(factor);
808 std::transform(std::begin(data), std::end(data), std::begin(data),
809 [coef](auto value) { return value * coef; });
810 }
811 }
812
813 void zero_centered_to_lb_in_place(auto &data) const {
815 }
816
817 void zero_centered_to_md_in_place(auto &data) const {
819 }
820
821 auto zero_centered_to_lb(auto const &data) const {
822 auto transformed_data = data;
823 zero_centered_to_lb_in_place(transformed_data);
824 return transformed_data;
825 }
826
827 auto zero_centered_to_md(auto const &data) const {
828 auto transformed_data = data;
829 zero_centered_to_md_in_place(transformed_data);
830 return transformed_data;
831 }
832
833public:
834 // ---- File I/O ----
835
836 void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override {
837 field::FlagFieldCellFilter<FlagField> fluid_filter(m_flag_field_id);
838 fluid_filter.addFlag(Boundary_flag);
839 vtk_obj.addCellExclusionFilter(fluid_filter);
840 }
841
842 void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj,
843 LatticeModel::units_map const &units,
844 int flag_observables) override;
845
846protected:
847 // ---- Private Infrastructure Helpers ----
848
849 /**
850 * @brief Convenience function to add a field with a custom allocator.
851 *
852 * When vectorization is off, let waLBerla decide which memory allocator
853 * to use. When vectorization is on, the aligned memory allocator is
854 * required, otherwise <tt>cpu_vectorize_info["assume_aligned"]</tt> will
855 * trigger assertions. That is because for single-precision kernels the
856 * waLBerla heuristic in <tt>src/field/allocation/FieldAllocator.h</tt>
857 * will fall back to @c StdFieldAlloc, yet @c AllocateAligned is needed
858 * for intrinsics to work.
859 */
860 template <typename Field> auto add_to_storage(std::string const tag) {
861 auto const &blocks = m_lattice->get_blocks();
862 auto const n_ghost_layers = m_lattice->get_ghost_layers();
863 if constexpr (Architecture == lbmpy::Arch::CPU) {
864#ifdef ESPRESSO_BUILD_WITH_AVX_KERNELS
865 constexpr auto alignment = field::SIMDAlignment();
866 using value_type = Field::value_type;
867 using Allocator = field::AllocateAligned<value_type, alignment>;
868 auto const allocator = std::make_shared<Allocator>();
869 auto const empty_set = Set<SUID>::emptySet();
870 return field::addToStorage<Field>(
871 blocks, tag, field::internal::defaultSize, FloatType{0}, field::fzyx,
872 n_ghost_layers, false, {}, empty_set, empty_set, allocator);
873#else // ESPRESSO_BUILD_WITH_AVX_KERNELS
874 return field::addToStorage<Field>(blocks, tag, FloatType{0}, field::fzyx,
875 n_ghost_layers);
876#endif // ESPRESSO_BUILD_WITH_AVX_KERNELS
877 }
878#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
879 else {
880 auto field_id = gpu::addGPUFieldToStorage<GPUField>(
881 blocks, tag, Field::F_SIZE, field::fzyx, n_ghost_layers);
882 if constexpr (std::is_same_v<Field, _VectorField>) {
883 for (auto &block : *blocks) {
884 auto field = block.template getData<GPUField>(field_id);
885 lbm::accessor::Vector::initialize(field, Vector3<FloatType>{0});
886 }
887 } else if constexpr (std::is_same_v<Field, _PdfField>) {
888 for (auto &block : *blocks) {
889 auto field = block.template getData<GPUField>(field_id);
891 field, std::array<FloatType, Stencil::Size>{});
892 }
893 }
894 return field_id;
895 }
896#endif
897 }
898
899 /**
900 * @brief Set up D3Q27 communicators for full ghost layer updates.
901 * Creates per-field communicators (PDF, velocity, last-applied force)
902 * as well as a combined communicator and the boundary communicator.
903 */
905 auto const &blocks = m_lattice->get_blocks();
906
907 m_full_communicator = std::make_shared<RegularFullCommunicator>(blocks);
908 m_full_communicator->addPackInfo(
909 std::make_shared<PackInfo<PdfField>>(m_pdf_field_id));
910 m_full_communicator->addPackInfo(
912 m_full_communicator->addPackInfo(
913 std::make_shared<PackInfo<VectorField>>(m_velocity_field_id));
914
915 m_pdf_communicator = std::make_shared<RegularFullCommunicator>(blocks);
916 m_vel_communicator = std::make_shared<RegularFullCommunicator>(blocks);
917 m_laf_communicator = std::make_shared<RegularFullCommunicator>(blocks);
918 m_pdf_communicator->addPackInfo(
919 std::make_shared<PackInfo<PdfField>>(m_pdf_field_id));
920 m_vel_communicator->addPackInfo(
921 std::make_shared<PackInfo<VectorField>>(m_velocity_field_id));
922 m_laf_communicator->addPackInfo(
924
926 std::make_shared<BoundaryFullCommunicator>(blocks);
927 m_boundary_communicator->addPackInfo(
930 auto boundary_packinfo = std::make_shared<
933 boundary_packinfo->setup_boundary_handle(m_lattice, m_boundary);
934 m_boundary_communicator->addPackInfo(boundary_packinfo);
935 }
936
937 /**
938 * @brief Set up the communicator used during integration.
939 * Uses optimized streaming pack info when neither boundaries nor
940 * Lees-Edwards boundary conditions are active; falls back to the
941 * generic pack info otherwise.
942 */
944 auto const setup = [this]<typename PackInfoPdf, typename PackInfoVec>() {
945 auto const &blocks = m_lattice->get_blocks();
947 std::make_shared<PDFStreamingCommunicator>(blocks);
948 m_pdf_streaming_communicator->addPackInfo(
949 std::make_shared<PackInfoPdf>(m_pdf_field_id));
950 m_pdf_streaming_communicator->addPackInfo(
951 std::make_shared<PackInfoVec>(m_last_applied_force_field_id));
952 };
954 using PackInfoPdf = FieldTrait::PackInfoStreamingPdf;
955 using PackInfoVec = FieldTrait::PackInfoStreamingVec;
956 if (m_has_boundaries or (m_collision_model and has_lees_edwards_bc())) {
957 setup.template operator()<PackInfo<PdfField>, PackInfoVec>();
958 } else {
959 setup.template operator()<PackInfoPdf, PackInfoVec>();
960 }
961 }
962};
963
964} // namespace walberla
965
966// Out-of-class template method definitions
970#include "LBNodeAccess.impl.hpp"
971#include "LBSliceAccess.impl.hpp"
972#include "LBVTK.impl.hpp"
Out-of-class boundary access definitions for walberla::LBWalberlaImpl.
Out-of-class collision model setup definitions for walberla::LBWalberlaImpl.
Out-of-class position-based interpolation definitions for walberla::LBWalberlaImpl.
Out-of-class node access definitions for walberla::LBWalberlaImpl.
Out-of-class slice access definitions for walberla::LBWalberlaImpl.
Out-of-class VTK writer registration definition for walberla::LBWalberlaImpl.
LBWalberlaBase provides the public interface of the LB waLBerla bridge.
Vector implementation and trait types for boost qvm interoperability.
Interface of a lattice-based fluid model.
std::map< std::string, std::shared_ptr< VTKHandle > > m_vtk_auto
VTK writers that are executed automatically.
std::unordered_map< std::string, double > units_map
Class that runs and controls the BlockForest in waLBerla.
walberla::blockforest::StructuredBlockForest Lattice_T
Boundary class optimized for sparse data.
field::FlagField< uint8_t > FlagField
Class that runs and controls the LB on waLBerla.
void add_forces_at_pos(std::vector< Utils::Vector3d > const &pos, std::vector< Utils::Vector3d > const &forces) override
Distribute forces to the lattice at given positions.
std::shared_ptr< typename Kernels::UpdateVelFromPDF > m_update_velocities_from_pdf
std::vector< double > get_slice_last_applied_force(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
std::variant< typename Kernels::StreamCollisionModelThermalized, typename Kernels::StreamCollisionModelLeesEdwards > CollisionModel
void zero_centered_transform_impl(T &data, auto const factor) const
Scale data by a conversion factor (in-place).
void reset_boundary_handling(std::shared_ptr< BlockStorage > const &blocks)
std::shared_ptr< RegularFullCommunicator > m_pdf_communicator
std::vector< std::optional< Utils::Vector3d > > get_slice_velocity_at_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
std::optional< Utils::Vector3d > get_node_last_applied_force(Utils::Vector3i const &node, bool consider_ghosts=false) const override
stencil::D3Q19 Stencil
Stencil for collision and streaming operations.
std::optional< Utils::Vector3d > get_node_velocity_at_boundary(Utils::Vector3i const &node, bool consider_ghosts=false) const override
void ghost_communication() override
Perform all pending ghost layer updates.
void pressure_tensor_correction(Matrix3< FloatType > &tensor) const
std::optional< Utils::Vector3d > get_node_velocity(Utils::Vector3i const &node, bool consider_ghosts=false) const override
FieldTrait< FloatType, Stencil, Architecture >::template RegularCommScheme< Stencil > PDFStreamingCommunicator
Regular communicator.
Utils::Vector3d get_boundary_force_from_shape(std::vector< int > const &raster_flat) const override
Total force exerted by the fluid on a subset of boundary nodes.
std::shared_ptr< RegularFullCommunicator > m_full_communicator
std::size_t get_force_field_id() const noexcept override
std::shared_ptr< CollisionModel > m_collision_model
void set_slice_velocity(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< double > const &velocity) override
unsigned int get_seed() const noexcept override
void integrate_vtk_writers() override
void pressure_tensor_correction(std::span< FloatType, 9ul > tensor) const
BoundaryModel::FlagField FlagField
bool remove_node_from_boundary(Utils::Vector3i const &node) override
std::optional< Utils::Vector3d > get_node_boundary_force(Utils::Vector3i const &node) const override
std::optional< double > get_density_at_pos(Utils::Vector3d const &pos, bool consider_points_in_halo=false) const override
void zero_centered_to_md_in_place(auto &data) const
FieldTrait< FloatType, Stencil >::VectorField _VectorField
void set_slice_velocity_at_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< std::optional< Utils::Vector3d > > const &velocity) override
std::vector< double > get_slice_population(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
void set_rng_state(uint64_t counter) override
std::vector< double > get_slice_velocity(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
Utils::Vector3d get_boundary_force() const override
std::shared_ptr< LatticeWalberla > m_lattice
std::shared_ptr< LeesEdwardsPack > m_lees_edwards_callbacks
LBWalberlaImpl(std::shared_ptr< LatticeWalberla > lattice, double viscosity, double density)
bool set_node_last_applied_force(Utils::Vector3i const &node, Utils::Vector3d const &force) override
std::function< bool(Utils::Vector3d const &)> make_lattice_position_checker(bool consider_points_in_halo) const override
FieldTrait< FloatType, Stencil >::PdfField _PdfField
std::shared_ptr< InterpolateAndShiftAtBoundary< _VectorField, FloatType > > m_lees_edwards_vel_interpol_sweep
std::vector< bool > get_slice_is_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
std::vector< double > get_slice_density(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
bool is_double_precision() const noexcept override
void update_boundary_from_shape(std::vector< int > const &raster_flat, std::vector< double > const &data_flat) override
Set boundary conditions from a rasterized shape.
std::shared_ptr< BoundaryModel > m_boundary
FieldTrait< FloatType, Stencil, Architecture >::PdfField PdfField
std::vector< double > get_slice_pressure_tensor(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
FieldTrait< FloatType, Stencil, Architecture >::template PackInfo< Field > PackInfo
bool set_node_density(Utils::Vector3i const &node, double density) override
void set_collision_model(double kT, unsigned int seed) override
Set up the thermalized collision model.
std::vector< Utils::Vector3d > get_velocities_at_pos(std::vector< Utils::Vector3d > const &pos) override
Interpolate velocities at given positions (batch version).
void ghost_communication_laf() override
std::optional< bool > get_node_is_boundary(Utils::Vector3i const &node, bool consider_ghosts=false) const override
std::shared_ptr< ResetForce< PdfField, VectorField > > m_reset_force
void setup_streaming_communicator()
Set up the communicator used during integration.
FlagUID const Boundary_flag
Flag for boundary cells.
~LBWalberlaImpl() override=default
std::optional< uint64_t > get_rng_state() const override
std::optional< std::vector< double > > get_node_population(Utils::Vector3i const &node, bool consider_ghosts=false) const override
FloatType m_density
kinematic viscosity
double get_viscosity() const noexcept override
double get_kT() const noexcept override
Utils::Vector3d get_momentum() const override
void set_viscosity(double viscosity) override
void on_boundary_add()
Lazily enable boundary mode on first boundary addition.
std::shared_ptr< BoundaryFullCommunicator > m_boundary_communicator
std::size_t stencil_size() const noexcept override
std::optional< BlockDataID > m_pressure_tensor_field_id
auto add_to_storage(std::string const tag)
Convenience function to add a field with a custom allocator.
void ghost_communication_full()
Communicate all fields at once using the D3Q27 stencil.
auto zero_centered_to_lb(auto const &data) const
FieldTrait< FloatType, Stencil, Architecture >::template BoundaryCommScheme< stencil::D3Q27 > BoundaryFullCommunicator
bool set_node_velocity_at_boundary(Utils::Vector3i const &node, Utils::Vector3d const &velocity) override
stencil::D3Q27 StencilFull
Stencil for ghost communication (includes domain corners).
void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj, LatticeModel::units_map const &units, int flag_observables) override
std::optional< double > get_node_density(Utils::Vector3i const &node, bool consider_ghosts=false) const override
void check_lebc(unsigned int shear_direction, unsigned int shear_plane_normal) const override
Verify that MD and LB Lees-Edwards parameters are consistent.
void set_slice_density(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< double > const &density) override
void ghost_communication_vel() override
Utils::Vector3d get_external_force() const noexcept override
FloatType FloatType_c(T t) const
std::shared_ptr< RegularFullCommunicator > m_laf_communicator
std::size_t get_velocity_field_id() const noexcept override
Utils::VectorXd< 9 > get_pressure_tensor() const override
ResourceObserver m_mpi_cart_comm_observer
FieldTrait< FloatType, Stencil, Architecture >::template RegularCommScheme< stencil::D3Q27 > RegularFullCommunicator
Full communicator.
std::shared_ptr< InterpolateAndShiftAtBoundary< _PdfField, FloatType > > m_lees_edwards_pdf_interpol_sweep
void set_external_force(Utils::Vector3d const &ext_force) override
void setup_full_communicator()
Set up D3Q27 communicators for full ghost layer updates.
std::optional< Utils::Vector3d > get_velocity_at_pos(Utils::Vector3d const &pos, bool consider_points_in_halo=false) const override
FieldTrait< FloatType, Stencil, Architecture >::VectorField VectorField
void set_slice_last_applied_force(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< double > const &force) override
LatticeWalberla::Lattice_T BlockStorage
Lattice model (e.g.
std::shared_ptr< RegularFullCommunicator > m_vel_communicator
std::bitset< GhostComm::SIZE > m_pending_ghost_comm
LatticeWalberla const & get_lattice() const noexcept override
auto zero_centered_to_md(auto const &data) const
void set_slice_population(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< double > const &population) override
double get_density() const noexcept override
void ghost_communication_pdf() override
void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override
std::optional< Utils::Vector3d > get_node_force_to_be_applied(Utils::Vector3i const &node) const override
std::vector< double > get_densities_at_pos(std::vector< Utils::Vector3d > const &pos) override
void zero_centered_to_lb_in_place(auto &data) const
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
bool is_gpu() const noexcept override
FloatType pressure_tensor_correction_factor() const
Correction factor for off-diagonal pressure tensor elements.
std::shared_ptr< PDFStreamingCommunicator > m_pdf_streaming_communicator
detail::KernelTrait< FloatType, Architecture > Kernels
BlockDataID m_last_applied_force_field_id
std::optional< Utils::VectorXd< 9 > > get_node_pressure_tensor(Utils::Vector3i const &node) const override
std::shared_ptr< InterpolateAndShiftAtBoundary< _VectorField, FloatType > > m_lees_edwards_last_applied_force_interpol_sweep
bool add_force_at_pos(Utils::Vector3d const &pos, Utils::Vector3d const &force) override
void setup_boundary_handle(std::shared_ptr< LatticeWalberla > lattice, std::shared_ptr< Boundary_T > boundary)
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:175
T product(Vector< T, N > const &v)
Definition Vector.hpp:372
STL namespace.
auto reduce(GhostLayerField< double, uint_t{19u}> const *pdf_field, GhostLayerField< double, uint_t{3u}> const *force_field, double const density)
void initialize(GhostLayerField< double, uint_t{19u}> *pdf_field, std::array< double, 19u > const &pop)
auto reduce(GhostLayerField< double, uint_t{19u}> const *pdf_field, double const density)
void initialize(GhostLayerField< double, uint_t{3u}> *vec_field, Vector3< double > const &vec)
\file PackInfoPdfDoublePrecision.cpp \author pystencils
auto to_vector3d(Vector3< T > const &v) noexcept
Cell to_cell(signed_integral_vector auto const &xyz)
ResourceObserver get_mpi_cart_comm_observer()
Get an observer on waLBerla's MPI Cartesian communicator status.
auto to_vector9d(Matrix3< T > const &m) noexcept
static Utils::Vector3d velocity(Particle const &p_ref, Particle const &p_vs)
Velocity of the virtual site.
Definition relative.cpp:64
Observer to monitor the lifetime of a shared resource.
detail::KernelTrait< FT, AT >::PackInfoVec PackInfoStreamingVec
Definition lb_fields.hpp:37
field::GhostLayerField< FT, PdfStencil::Size > PdfField
Definition lb_fields.hpp:33
field::GhostLayerField< FT, uint_t{3u}> VectorField
Definition lb_fields.hpp:34
detail::KernelTrait< FT, AT >::PackInfoPdf PackInfoStreamingPdf
Definition lb_fields.hpp:36
GhostCommFlags
Ghost communication operations.
@ LAF
last applied forces communication