ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
EKinWalberlaImpl.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2022-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#include <blockforest/communication/UniformBufferedScheme.h>
23#include <field/AddToStorage.h>
24#include <field/FlagField.h>
25#include <field/FlagUID.h>
26#include <field/GhostLayerField.h>
27#include <field/communication/PackInfo.h>
28#include <field/iterators/IteratorMacros.h>
29#include <field/vtk/FlagFieldCellFilter.h>
30#include <field/vtk/VTKWriter.h>
31#include <stencil/D3Q27.h>
32#include <waLBerlaDefinitions.h>
33#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
34#include <gpu/AddGPUFieldToStorage.h>
35#include <gpu/communication/MemcpyPackInfo.h>
36#include <gpu/communication/UniformGPUScheme.h>
37#endif
38
39#include "../BoundaryHandling.hpp"
40#include "../BoundaryPackInfo.hpp"
41#include "../utils/boundary.hpp"
42#include "../utils/types_conversion.hpp"
43#include "ek_kernels.hpp"
44#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
45#include "ek_kernels.cuh"
46#endif
47
54
55#include <utils/Vector.hpp>
56
57#include <cstddef>
58#include <iterator>
59#include <memory>
60#include <optional>
61#include <ranges>
62#include <stdexcept>
63#include <string>
64#include <type_traits>
65#include <variant>
66#include <vector>
67
68namespace walberla {
69
70/** @brief Class that runs and controls the EK on waLBerla. */
71template <std::size_t FluxCount = 13, typename FloatType = double,
74#if not defined(WALBERLA_BUILD_WITH_CUDA)
75 static_assert(Architecture != lbmpy::Arch::GPU,
76 "waLBerla was compiled without CUDA support");
77#endif
78 using ContinuityKernel =
80 using DiffusiveFluxKernelUnthermalized =
81 typename detail::KernelTrait<FloatType,
82 Architecture>::DiffusiveFluxKernel;
83 using DiffusiveFluxKernelThermalized = typename detail::KernelTrait<
84 FloatType, Architecture>::DiffusiveFluxKernelThermalized;
85 using AdvectiveFluxKernel =
86 typename detail::KernelTrait<FloatType,
87 Architecture>::AdvectiveFluxKernel;
88 using FrictionCouplingKernel =
89 typename detail::KernelTrait<FloatType,
90 Architecture>::FrictionCouplingKernel;
91 using DiffusiveFluxKernelElectrostaticUnthermalized =
92 typename detail::KernelTrait<
93 FloatType, Architecture>::DiffusiveFluxKernelElectrostatic;
94 using DiffusiveFluxKernelElectrostaticThermalized =
95 typename detail::KernelTrait<
96 FloatType, Architecture>::DiffusiveFluxKernelElectrostaticThermalized;
97
98 using DiffusiveFluxKernel = std::variant<DiffusiveFluxKernelUnthermalized,
99 DiffusiveFluxKernelThermalized>;
100 using DiffusiveFluxKernelElectrostatic =
101 std::variant<DiffusiveFluxKernelElectrostaticUnthermalized,
102 DiffusiveFluxKernelElectrostaticThermalized>;
103
104 using Dirichlet =
106 using FixedFlux =
108
111 using BoundaryModelFlux =
113
114public:
115 /** @brief Stencil for collision and streaming operations. */
116 using Stencil = stencil::D3Q27;
117 /** @brief Lattice model (e.g. blockforest). */
119
120protected:
121 template <typename FT, lbmpy::Arch AT = lbmpy::Arch::CPU> struct FieldTrait {
122 // Type definitions
125 template <class Field>
126 using PackInfo = field::communication::PackInfo<Field>;
127 template <class Stencil>
129 blockforest::communication::UniformBufferedScheme<Stencil>;
130 template <class Stencil>
132 blockforest::communication::UniformBufferedScheme<Stencil>;
133 };
134 using FlagField = walberla::FlagField<walberla::uint8_t>;
135#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
136 template <typename FT> struct FieldTrait<FT, lbmpy::Arch::GPU> {
137 private:
138 static auto constexpr AT = lbmpy::Arch::GPU;
139 template <class Field>
140 using MemcpyPackInfo = gpu::communication::MemcpyPackInfo<Field>;
141
142 public:
143 template <typename Stencil>
144 class UniformGPUScheme
145 : public gpu::communication::UniformGPUScheme<Stencil> {
146 public:
147 explicit UniformGPUScheme(auto const &bf)
148 : gpu::communication::UniformGPUScheme<Stencil>(
149 bf, /* sendDirectlyFromGPU */ false,
150 /* useLocalCommunication */ false) {}
151 };
152 using FluxField = gpu::GPUField<FT>;
153 using DensityField = gpu::GPUField<FT>;
154 template <class Field> using PackInfo = MemcpyPackInfo<Field>;
155 template <class Stencil>
157 template <class Stencil>
158 using BoundaryCommScheme =
159 blockforest::communication::UniformBufferedScheme<Stencil>;
160 };
161 using GPUField = gpu::GPUField<FloatType>;
162#endif
163
164 struct GhostComm {
165 /** @brief Ghost communication operations. */
166 enum GhostCommFlags : unsigned {
167 FLB, ///< flux boundary communication
168 DENS, ///< density communication
169 SIZE
170 };
171 };
172
173 // "underlying" field types (`GPUField` has no f-size info at compile time)
176
177public:
181
182 template <typename T> FloatType FloatType_c(T t) {
183 return numeric_cast<FloatType>(t);
184 }
185
186 [[nodiscard]] std::size_t stencil_size() const noexcept override {
187 return FluxCount;
188 }
189
191 return std::is_same_v<FloatType, double>;
192 }
193
194private:
195 FloatType m_diffusion;
196 FloatType m_kT;
197 FloatType m_valency;
198 Utils::Vector3d m_ext_efield;
199 bool m_advection;
200 bool m_friction_coupling;
201 unsigned int m_seed;
202
203protected:
204 // Block data access handles
206
208
211
212 /** Flag for domain cells, i.e. all cells. */
213 FlagUID const Domain_flag{"domain"};
214 /** Flag for boundary cells. */
215 FlagUID const Boundary_flag{"boundary"};
216
217 /** Block forest */
218 std::shared_ptr<LatticeWalberla> m_lattice;
219
220 std::unique_ptr<BoundaryModelDensity> m_boundary_density;
221 std::shared_ptr<BoundaryModelFlux> m_boundary_flux;
222
223 std::unique_ptr<DiffusiveFluxKernel> m_diffusive_flux;
224 std::unique_ptr<DiffusiveFluxKernelElectrostatic>
226 std::unique_ptr<ContinuityKernel> m_continuity;
227
228 // ResetFlux + external force
229 // TODO: kernel for that
230 // std::shared_ptr<ResetForce<PdfField, VectorField>> m_reset_force;
231
232 /**
233 * @brief Convenience function to add a field with a custom allocator.
234 *
235 * When vectorization is off, let waLBerla decide which memory allocator
236 * to use. When vectorization is on, the aligned memory allocator is
237 * required, otherwise <tt>cpu_vectorize_info["assume_aligned"]</tt> will
238 * trigger assertions. That is because for single-precision kernels the
239 * waLBerla heuristic in <tt>src/field/allocation/FieldAllocator.h</tt>
240 * will fall back to @c StdFieldAlloc, yet @c AllocateAligned is needed
241 * for intrinsics to work.
242 */
243 template <typename Field>
244 auto add_to_storage(std::string const tag, FloatType value) {
245 auto const &blocks = m_lattice->get_blocks();
246 auto const n_ghost_layers = m_lattice->get_ghost_layers();
247#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
248 if constexpr (Architecture == lbmpy::Arch::GPU) {
249 auto field_id = gpu::addGPUFieldToStorage<GPUField>(
250 blocks, tag, Field::F_SIZE, field::fzyx, n_ghost_layers);
251 if constexpr (std::is_same_v<Field, _DensityField>) {
252 for (auto block = blocks->begin(); block != blocks->end(); ++block) {
253 auto field = block->template getData<GPUField>(field_id);
254 ek::accessor::Scalar::initialize(field, FloatType{value});
255 }
256 } else if constexpr (std::is_same_v<Field, _FluxField>) {
257 for (auto block = blocks->begin(); block != blocks->end(); ++block) {
258 auto field = block->template getData<GPUField>(field_id);
260 std::array<FloatType, FluxCount>{});
261 }
262 }
263 return field_id;
264 }
265#endif
266 return field::addToStorage<Field>(blocks, tag, FloatType{value},
267 field::fzyx, n_ghost_layers);
268 }
269
270 void
271 reset_density_boundary_handling(std::shared_ptr<BlockStorage> const &blocks) {
272 auto const [lc, uc] = m_lattice->get_local_grid_range(true);
273 m_boundary_density = std::make_unique<BoundaryModelDensity>(
275 CellInterval{to_cell(lc), to_cell(uc)});
276 }
277
278 void
279 reset_flux_boundary_handling(std::shared_ptr<BlockStorage> const &blocks) {
280 auto const [lc, uc] = m_lattice->get_local_grid_range(true);
281 m_boundary_flux = std::make_shared<BoundaryModelFlux>(
283 CellInterval{to_cell(lc), to_cell(uc)});
284 }
285
287 typename FieldTrait<FloatType, Architecture>::template RegularCommScheme<
288 typename stencil::D3Q27>;
290 typename FieldTrait<FloatType, Architecture>::template BoundaryCommScheme<
291 typename stencil::D3Q27>;
292 std::shared_ptr<FullCommunicator> m_full_communication;
293 std::shared_ptr<BoundaryFullCommunicator> m_boundary_communicator;
294 std::bitset<GhostComm::SIZE> m_pending_ghost_comm;
296 template <class Field>
297 using PackInfo =
299
300public:
301 EKinWalberlaImpl(std::shared_ptr<LatticeWalberla> lattice, double diffusion,
302 double kT, double valency, Utils::Vector3d const &ext_efield,
303 double density, bool advection, bool friction_coupling,
304 bool thermalized, unsigned int seed)
305 : m_diffusion(FloatType_c(diffusion)), m_kT(FloatType_c(kT)),
306 m_valency(FloatType_c(valency)), m_ext_efield(ext_efield),
307 m_advection(advection), m_friction_coupling(friction_coupling),
308 m_seed(seed), m_lattice(std::move(lattice)),
310
311 auto const &blocks = m_lattice->get_blocks();
312 auto const n_ghost_layers = m_lattice->get_ghost_layers();
313
317 add_to_storage<_FluxField>("flux field", FloatType_c(0.0));
318
320 std::make_unique<ContinuityKernel>(m_flux_field_id, m_density_field_id);
321
322 if (thermalized) {
323 set_diffusion_kernels(*m_lattice, seed);
324 } else {
325 set_diffusion_kernels();
326 }
327
328 // Init boundary related stuff
329 m_flag_field_density_id = field::addFlagFieldToStorage<FlagField>(
330 blocks, "flag field density", n_ghost_layers);
332
333 m_flag_field_flux_id = field::addFlagFieldToStorage<FlagField>(
334 blocks, "flag field flux", n_ghost_layers);
336
337 m_full_communication = std::make_shared<FullCommunicator>(blocks);
338 m_full_communication->addPackInfo(
339 std::make_shared<PackInfo<DensityField>>(m_density_field_id));
341 std::make_shared<BoundaryFullCommunicator>(blocks);
342 m_boundary_communicator->addPackInfo(
345 auto flux_boundary_packinfo = std::make_shared<
350
352 }
353
354 // Global parameters
355 [[nodiscard]] double get_diffusion() const noexcept override {
356 return m_diffusion;
357 }
358 [[nodiscard]] double get_kT() const noexcept override { return m_kT; }
359 [[nodiscard]] double get_valency() const noexcept override {
360 return m_valency;
361 }
362 [[nodiscard]] bool get_advection() const noexcept override {
363 return m_advection;
364 }
366 return m_friction_coupling;
367 }
369 return m_ext_efield;
370 }
372 return static_cast<bool>(
373 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux));
374 }
375 [[nodiscard]] unsigned int get_seed() const noexcept override {
376 return m_seed;
377 }
378 [[nodiscard]] std::optional<uint64_t> get_rng_state() const override {
379 auto const kernel =
380 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux);
381 if (!kernel) {
382 return std::nullopt;
383 }
384 return {static_cast<uint64_t>(kernel->getTime_step())};
385 }
386
387 void set_diffusion(double diffusion) override {
388 m_diffusion = FloatType_c(diffusion);
389 auto visitor = [m_diffusion = m_diffusion](auto &kernel) {
390 kernel.setD(m_diffusion);
391 };
392 std::visit(visitor, *m_diffusive_flux);
394 }
395
396 void set_kT(double kT) override {
397 m_kT = FloatType_c(kT);
398 std::visit([m_kT = m_kT](auto &kernel) { kernel.setKt(m_kT); },
400 }
401
402 void set_valency(double valency) override {
403 m_valency = FloatType_c(valency);
404 std::visit(
405 [m_valency = m_valency](auto &kernel) { kernel.setZ(m_valency); },
407 }
408
409 void set_advection(bool advection) override { m_advection = advection; }
410
412 m_friction_coupling = friction_coupling;
413 }
414
415 void set_rng_state(uint64_t counter) override {
416 auto const kernel =
417 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux);
418 auto const kernel_electrostatic =
419 std::get_if<DiffusiveFluxKernelElectrostaticThermalized>(
421
422 if (!kernel or !kernel_electrostatic) {
423 throw std::runtime_error("This EK instance is unthermalized");
424 }
425 assert(counter <=
426 static_cast<uint32_t>(std::numeric_limits<uint_t>::max()));
427 kernel->setTime_step(static_cast<uint32_t>(counter));
428 kernel_electrostatic->setTime_step(static_cast<uint32_t>(counter));
429 }
430
431 void set_ext_efield(Utils::Vector3d const &field) override {
432 m_ext_efield = field;
433
434 std::visit(
435 [this](auto &kernel) {
436 kernel.setF_ext_0(FloatType_c(m_ext_efield[0]));
437 kernel.setF_ext_1(FloatType_c(m_ext_efield[1]));
438 kernel.setF_ext_2(FloatType_c(m_ext_efield[2]));
439 },
441 }
442
443 void ghost_communication() override {
446 (*m_full_communication)();
448 }
450 }
451
459
460private:
461 void set_diffusion_kernels() {
462 auto kernel = DiffusiveFluxKernelUnthermalized(
464 m_diffusive_flux = std::make_unique<DiffusiveFluxKernel>(std::move(kernel));
465
466 auto kernel_electrostatic = DiffusiveFluxKernelElectrostaticUnthermalized(
468 FloatType_c(m_diffusion), FloatType_c(m_ext_efield[0]),
469 FloatType_c(m_ext_efield[1]), FloatType_c(m_ext_efield[2]),
470 FloatType_c(m_kT), FloatType_c(m_valency));
471
473 std::make_unique<DiffusiveFluxKernelElectrostatic>(
474 std::move(kernel_electrostatic));
475 }
476
477 void set_diffusion_kernels(LatticeWalberla const &lattice,
478 unsigned int seed) {
479 auto const grid_dim = lattice.get_grid_dimensions();
480
481 auto kernel = DiffusiveFluxKernelThermalized(
483 grid_dim[0], grid_dim[1], grid_dim[2], seed, 0);
484
485 auto kernel_electrostatic = DiffusiveFluxKernelElectrostaticThermalized(
487 FloatType_c(m_diffusion), FloatType_c(m_ext_efield[0]),
488 FloatType_c(m_ext_efield[1]), FloatType_c(m_ext_efield[2]), grid_dim[0],
489 grid_dim[1], grid_dim[2], FloatType_c(m_kT), seed, 0,
490 FloatType_c(m_valency));
491
492 auto const blocks = lattice.get_blocks();
493
494 for (auto &block : *blocks) {
495 kernel.configure(blocks, &block);
496 kernel_electrostatic.configure(blocks, &block);
497 }
498
499 m_diffusive_flux = std::make_unique<DiffusiveFluxKernel>(std::move(kernel));
501 std::make_unique<DiffusiveFluxKernelElectrostatic>(
502 std::move(kernel_electrostatic));
503 }
504
505 void kernel_boundary_density() {
506 for (auto &block : *m_lattice->get_blocks()) {
507 (*m_boundary_density)(&block);
508 }
509 }
510
511 void kernel_boundary_flux() {
512 for (auto &block : *m_lattice->get_blocks()) {
513 (*m_boundary_flux)(&block);
514 }
515 }
516
517 void kernel_continuity() {
518 for (auto &block : *m_lattice->get_blocks()) {
519 (*m_continuity).run(&block);
520 }
521 }
522
523 void kernel_diffusion() {
524 for (auto &block : *m_lattice->get_blocks()) {
525 std::visit([&block](auto &kernel) { kernel.run(&block); },
527 }
528
529 if (auto *kernel =
530 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux)) {
531 kernel->setTime_step(kernel->getTime_step() + 1u);
532
534 std::get_if<DiffusiveFluxKernelElectrostaticThermalized>(
536 kernel_electrostatic->setTime_step(kernel_electrostatic->getTime_step() +
537 1u);
538 }
539 }
540
541 void kernel_advection(std::size_t const velocity_id) {
542 auto kernel = AdvectiveFluxKernel(m_flux_field_id, m_density_field_id,
544 for (auto &block : *m_lattice->get_blocks()) {
545 kernel.run(&block);
546 }
547 }
548
549 void kernel_friction_coupling(std::size_t const force_id,
550 double const lb_density) {
551 auto kernel = FrictionCouplingKernel(
553 FloatType_c(get_kT()), FloatType(lb_density));
554 for (auto &block : *m_lattice->get_blocks()) {
555 kernel.run(&block);
556 }
557 }
558
559 void kernel_diffusion_electrostatic(std::size_t const potential_id) {
560 auto const phiID = BlockDataID(potential_id);
561 std::visit([phiID](auto &kernel) { kernel.setPhiID(phiID); },
563
564 for (auto &block : *m_lattice->get_blocks()) {
565 std::visit([&block](auto &kernel) { kernel.run(&block); },
567 }
568
569 if (auto *kernel_electrostatic =
570 std::get_if<DiffusiveFluxKernelElectrostaticThermalized>(
572 kernel_electrostatic->setTime_step(kernel_electrostatic->getTime_step() +
573 1u);
574
575 auto *kernel =
576 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux);
577 kernel->setTime_step(kernel->getTime_step() + 1u);
578 }
579 }
580
581 void kernel_migration() {}
582
583 void update_boundary_fields() {
584 m_boundary_flux->boundary_update();
585 m_boundary_density->boundary_update();
586 }
587
588protected:
589 void integrate_vtk_writers() override {
590 for (auto const &vtk_handle : m_vtk_auto | std::views::values) {
591 if (vtk_handle->enabled) {
592 vtk::writeFiles(vtk_handle->ptr)();
593 vtk_handle->execution_count++;
594 }
595 }
596 }
597
598public:
599 void integrate(std::size_t potential_id, std::size_t velocity_id,
600 std::size_t force_id, double lb_density) override {
601
603 update_boundary_fields();
604
605 if (get_diffusion() == 0.)
606 return;
607
608 if (get_valency() != 0.) {
609 if (potential_id == walberla::BlockDataID{}) {
610 throw std::runtime_error("Walberla EK: electrostatic potential enabled "
611 "but no field accessible. potential id is " +
612 std::to_string(potential_id));
613 }
614 kernel_diffusion_electrostatic(potential_id);
615 } else {
616 kernel_diffusion();
617 }
618
619 kernel_migration();
620 kernel_boundary_flux();
621 // friction coupling
622 if (get_friction_coupling()) {
623 if (force_id == walberla::BlockDataID{}) {
624 throw std::runtime_error("Walberla EK: friction coupling enabled but "
625 "no force field accessible. force_id is " +
626 std::to_string(force_id) +
627 ". Hint: LB may be inactive.");
628 }
629 kernel_friction_coupling(force_id, lb_density);
630 }
631
632 if (get_advection()) {
633 if (velocity_id == walberla::BlockDataID{}) {
634 throw std::runtime_error("Walberla EK: advection enabled but no "
635 "velocity field accessible. velocity_id is " +
636 std::to_string(velocity_id) +
637 ". Hint: LB may be inactive.");
638 }
639 kernel_advection(velocity_id);
640 kernel_boundary_flux();
641 }
642 kernel_continuity();
643
644 // is this the expected behavior when reactions are included?
645 kernel_boundary_density();
647
648 // Handle VTK writers
650 }
651
652 [[nodiscard]] std::size_t get_density_id() const noexcept override {
653 static_assert(std::is_same_v<std::size_t, walberla::uint_t>);
654 return static_cast<std::size_t>(m_density_field_id);
655 }
656
657 bool set_node_density(Utils::Vector3i const &node, double density) override {
659 auto bc = get_block_and_cell(get_lattice(), node, false);
660 if (!bc)
661 return false;
662
663 auto density_field =
666 return true;
667 }
668
669 [[nodiscard]] std::optional<double>
671 bool consider_ghosts = false) const override {
672 if (m_boundary_density->node_is_boundary(node)) {
673 return m_boundary_density->get_node_value_at_boundary(node);
674 }
675
677
678 if (!bc)
679 return std::nullopt;
680
681 auto const density_field =
684 }
685
686 [[nodiscard]] std::vector<double>
688 Utils::Vector3i const &upper_corner) const override {
689 std::vector<double> out;
690#ifndef NDEBUG
692#endif
693 auto const &lattice = get_lattice();
694 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
695 out = std::vector<double>(ci->numCells());
696 for (auto &block : *lattice.get_blocks()) {
697 auto const block_offset = lattice.get_block_corner(block, true);
698 if (auto const bci = get_block_interval(
700 auto const density_field =
703 assert(values.size() == bci->numCells());
704#ifndef NDEBUG
705 values_size += bci->numCells();
706#endif
707 auto kernel = [this, &values, &out](unsigned const block_index,
708 unsigned const local_index,
709 Utils::Vector3i const &node) {
710 if (m_boundary_density->node_is_boundary(node)) {
712 m_boundary_density->get_node_value_at_boundary(node);
713 } else {
715 }
716 };
717
719 }
720 }
721 assert(values_size == ci->numCells());
722 }
723 return out;
724 }
725
728 std::vector<double> const &density) override {
730 auto const &lattice = get_lattice();
731 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
732 assert(density.size() == ci->numCells());
733 for (auto &block : *lattice.get_blocks()) {
734 auto const block_offset = lattice.get_block_corner(block, true);
735 if (auto const bci = get_block_interval(
737 auto const density_field =
739 std::vector<FloatType> values(bci->numCells());
740
741 auto kernel = [&values, &density](unsigned const block_index,
742 unsigned const local_index,
743 Utils::Vector3i const &) {
745 };
746
749 }
750 }
751 }
752 }
753
754 [[nodiscard]] std::optional<Utils::Vector3d>
756 bool consider_ghosts = false) const override {
757 if (m_boundary_flux->node_is_boundary(node)) {
758 return to_vector3d(m_boundary_flux->get_node_value_at_boundary(node));
759 }
760
762
763 if (!bc)
764 return std::nullopt;
765
766 auto const flux_field =
767 bc->block->template getData<FluxField>(m_flux_field_id);
769 }
770
771 std::vector<double>
773 Utils::Vector3i const &upper_corner) const override {
774 std::vector<double> out;
775#ifndef NDEBUG
777#endif
778 auto const &lattice = get_lattice();
779 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
780 out = std::vector<double>(3u * ci->numCells());
781 for (auto &block : *lattice.get_blocks()) {
782 auto const block_offset = lattice.get_block_corner(block, true);
783 if (auto const bci = get_block_interval(
785 auto const flux_field =
788 assert(values.size() == 3u * bci->numCells());
789#ifndef NDEBUG
790 values_size += 3u * bci->numCells();
791#endif
792
793 auto kernel = [&values, &out, this](unsigned const block_index,
794 unsigned const local_index,
795 Utils::Vector3i const &node) {
796 if (m_boundary_flux->node_is_boundary(node)) {
797 auto const &vec =
798 m_boundary_flux->get_node_value_at_boundary(node);
799 for (uint_t f = 0u; f < 3u; ++f) {
800 out[3u * local_index + f] = double_c(vec[f]);
801 }
802 } else {
803 for (uint_t f = 0u; f < 3u; ++f) {
804 out[3u * local_index + f] =
805 double_c(values[3u * block_index + f]);
806 }
807 }
808 };
809
811 }
812 }
813 assert(values_size == 3u * ci->numCells());
814 }
815 return out;
816 }
817
822
823 void clear_density_boundaries() override {
825 }
826
828 Utils::Vector3d const &flux) override {
830 auto bc = get_block_and_cell(get_lattice(), node, true);
831 if (!bc)
832 return false;
833
834 m_boundary_flux->set_node_value_at_boundary(
835 node, to_vector3<FloatType>(flux), *bc);
836 return true;
837 }
838
839 [[nodiscard]] std::optional<Utils::Vector3d>
841 bool consider_ghosts = false) const override {
843 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
844 if (!bc or !m_boundary_flux->node_is_boundary(node))
845 return std::nullopt;
846
847 return {to_vector3d(m_boundary_flux->get_node_value_at_boundary(node))};
848 }
849
852 auto bc = get_block_and_cell(get_lattice(), node, true);
853 if (!bc)
854 return false;
855
856 m_boundary_flux->remove_node_from_boundary(node, *bc);
857 return true;
858 }
859
861 double density) override {
862 auto bc = get_block_and_cell(get_lattice(), node, true);
863 if (!bc)
864 return false;
865
866 m_boundary_density->set_node_value_at_boundary(node, FloatType_c(density),
867 *bc);
868
869 return true;
870 }
871
872 [[nodiscard]] std::optional<double>
874 bool consider_ghosts = false) const override {
875 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
876 if (!bc or !m_boundary_density->node_is_boundary(node))
877 return std::nullopt;
878
879 return {double_c(m_boundary_density->get_node_value_at_boundary(node))};
880 }
881
884 std::vector<std::optional<double>> const &density) override {
885 auto const &lattice = get_lattice();
886 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
887 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
888 auto const lower_cell = ci->min();
889 auto const upper_cell = ci->max();
890 auto it = density.begin();
891 assert(density.size() == ci->numCells());
892 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
893 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
894 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
895 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
896 auto const bc = get_block_and_cell(lattice, node, false);
897 auto const &opt = *it;
898 if (opt) {
899 m_boundary_density->set_node_value_at_boundary(
900 node, FloatType_c(*opt), *bc);
901 } else {
902 m_boundary_density->remove_node_from_boundary(node, *bc);
903 }
904 ++it;
905 }
906 }
907 }
908 }
909 }
910
911 [[nodiscard]] std::vector<std::optional<double>>
914 Utils::Vector3i const &upper_corner) const override {
915 std::vector<std::optional<double>> out;
916 auto const &lattice = get_lattice();
917 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
918 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
919 auto const lower_cell = ci->min();
920 auto const upper_cell = ci->max();
921 auto const n_values = ci->numCells();
922 out.reserve(n_values);
923 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
924 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
925 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
926 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
927 if (m_boundary_density->node_is_boundary(node)) {
928 out.emplace_back(double_c(
929 m_boundary_density->get_node_value_at_boundary(node)));
930 } else {
931 out.emplace_back(std::nullopt);
932 }
933 }
934 }
935 }
936 assert(out.size() == n_values);
937 }
938 return out;
939 }
940
943 std::vector<std::optional<Utils::Vector3d>> const &flux) override {
945 auto const &lattice = get_lattice();
946 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
947 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
948 auto const lower_cell = ci->min();
949 auto const upper_cell = ci->max();
950 auto it = flux.begin();
951 assert(flux.size() == ci->numCells());
952 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
953 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
954 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
955 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
956 auto const bc = get_block_and_cell(lattice, node, false);
957 auto const &opt = *it;
958 if (opt) {
959 m_boundary_flux->set_node_value_at_boundary(
960 node, to_vector3<FloatType>(*opt), *bc);
961 } else {
962 m_boundary_flux->remove_node_from_boundary(node, *bc);
963 }
964 ++it;
965 }
966 }
967 }
968 }
969 }
970
971 [[nodiscard]] std::vector<std::optional<Utils::Vector3d>>
974 Utils::Vector3i const &upper_corner) const override {
975 std::vector<std::optional<Utils::Vector3d>> out;
976 auto const &lattice = get_lattice();
977 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
978 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
979 auto const lower_cell = ci->min();
980 auto const upper_cell = ci->max();
981 auto const n_values = ci->numCells();
982 out.reserve(n_values);
983 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
984 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
985 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
986 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
987 if (m_boundary_flux->node_is_boundary(node)) {
988 out.emplace_back(to_vector3d(
989 m_boundary_flux->get_node_value_at_boundary(node)));
990 } else {
991 out.emplace_back(std::nullopt);
992 }
993 }
994 }
995 }
996 assert(out.size() == n_values);
997 }
998 return out;
999 }
1000
1001 [[nodiscard]] std::vector<bool>
1003 Utils::Vector3i const &upper_corner) const override {
1004 std::vector<bool> out;
1005 auto const &lattice = get_lattice();
1006 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
1007 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
1008 auto const lower_cell = ci->min();
1009 auto const upper_cell = ci->max();
1010 auto const n_values = ci->numCells();
1011 out.reserve(n_values);
1012 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
1013 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
1014 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
1015 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
1016 out.emplace_back(m_boundary_density->node_is_boundary(node) or
1017 m_boundary_flux->node_is_boundary(node));
1018 }
1019 }
1020 }
1021 assert(out.size() == n_values);
1022 }
1023 return out;
1024 }
1025
1027 auto bc = get_block_and_cell(get_lattice(), node, true);
1028 if (!bc)
1029 return false;
1030
1031 m_boundary_density->remove_node_from_boundary(node, *bc);
1032
1033 return true;
1034 }
1035
1036 [[nodiscard]] std::optional<bool>
1038 bool consider_ghosts) const override {
1041 if (!bc)
1042 return std::nullopt;
1043
1044 return {m_boundary_flux->node_is_boundary(node)};
1045 }
1046
1047 [[nodiscard]] std::optional<bool>
1049 bool consider_ghosts) const override {
1051 if (!bc)
1052 return std::nullopt;
1053
1054 return {m_boundary_density->node_is_boundary(node)};
1055 }
1056
1057 [[nodiscard]] std::optional<bool>
1059 bool consider_ghosts = false) const override {
1061 if (!bc)
1062 return std::nullopt;
1063
1064 return {m_boundary_density->node_is_boundary(node) or
1065 m_boundary_flux->node_is_boundary(node)};
1066 }
1067
1069 const std::vector<int> &raster_flat,
1070 const std::vector<double> &data_flat) override {
1072 auto const grid_size = get_lattice().get_grid_dimensions();
1073 auto const data = fill_3D_vector_array(data_flat, grid_size);
1076 }
1077
1079 const std::vector<int> &raster_flat,
1080 const std::vector<double> &data_flat) override {
1081 auto const grid_size = get_lattice().get_grid_dimensions();
1082 auto const data = fill_3D_scalar_array(data_flat, grid_size);
1084 data);
1086 }
1087
1089
1091 m_boundary_density->boundary_update();
1092 }
1093
1095 return *m_lattice;
1096 }
1097
1098 [[nodiscard]] bool is_gpu() const noexcept override {
1100 }
1101
1102 void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override {
1103 field::FlagFieldCellFilter<FlagField> dens_filter(m_flag_field_density_id);
1104 field::FlagFieldCellFilter<FlagField> flux_filter(m_flag_field_flux_id);
1105 dens_filter.addFlag(Boundary_flag);
1106 flux_filter.addFlag(Boundary_flag);
1107 vtk_obj.addCellExclusionFilter(dens_filter);
1108 vtk_obj.addCellExclusionFilter(flux_filter);
1109 }
1110
1111protected:
1112 template <typename VecType, uint_t F_SIZE_ARG, typename OutputType>
1113 class VTKWriter : public vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG> {
1114 public:
1115 VTKWriter(ConstBlockDataID const &block_id, std::string const &id,
1116 FloatType unit_conversion)
1117 : vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG>(id),
1119
1120 protected:
1122
1123 std::size_t get_first_index(cell_idx_t const x, cell_idx_t const y,
1124 cell_idx_t const z) {
1125 return (static_cast<std::size_t>(x) * m_dims[2] * m_dims[1] +
1126 static_cast<std::size_t>(y) * m_dims[2] +
1127 static_cast<std::size_t>(z)) *
1128 F_SIZE_ARG;
1129 }
1130
1131 FloatType m_conversion;
1134
1135 public:
1137
1138 void set_dims(Vector3<uint_t> dims) { m_dims = dims; }
1139 };
1140
1141 template <typename OutputType = float>
1143 : public VTKWriter<std::vector<FloatType>, 1u, OutputType> {
1144 public:
1146 using Base::Base;
1147 using Base::evaluate;
1148
1149 protected:
1151 cell_idx_t const z, cell_idx_t const) override {
1152 WALBERLA_ASSERT(!this->m_content.empty());
1153 auto const density = this->m_content[this->get_first_index(x, y, z)];
1155 }
1156 };
1157
1158 template <typename OutputType = float>
1160 : public VTKWriter<std::vector<FloatType>, 3u, OutputType> {
1161 public:
1163 using Base::Base;
1164 using Base::evaluate;
1165
1166 protected:
1168 cell_idx_t const z, cell_idx_t const f) override {
1169 WALBERLA_ASSERT(!this->m_content.empty());
1170 auto flux = this->m_content[this->get_first_index(x, y, z) + f];
1172 }
1173 };
1174
1175 template <typename OutputType = float>
1176 class BoundaryVTKWriter : public vtk::BlockCellDataWriter<OutputType, 1u> {
1177 public:
1178 using Base = vtk::BlockCellDataWriter<OutputType, 1u>;
1179 using Base::evaluate;
1181 std::string const &id, FlagUID const &boundary_flag)
1182 : vtk::BlockCellDataWriter<OutputType, 1u>(id),
1185
1186 protected:
1192
1194 cell_idx_t const z, cell_idx_t const) override {
1196 return m_flag_field->isFlagSet(x, y, z, m_boundary_flag_value)
1197 ? OutputType{1}
1198 : OutputType{0};
1199 }
1200
1204 typename FlagField::flag_t m_boundary_flag_value;
1205 };
1206
1207public:
1208 void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj,
1210 int flag_observables) override {
1211 if (flag_observables & static_cast<int>(EKOutputVTK::density)) {
1212 auto const unit_conversion = FloatType_c(units.at("density"));
1213 auto const blocks = m_lattice->get_blocks();
1217 auto before_function = [this, blocks, density_writer]() {
1218 for (auto &block : *blocks) {
1219 auto *density_field =
1221
1222 auto const offset = m_lattice->get_block_corner(block, true);
1223 auto const *flag_field =
1225 auto const boundary_flag = flag_field->getFlag(Boundary_flag);
1227 if (flag_field->isFlagSet(x, y, z, boundary_flag)) {
1228 Cell const global(offset[0] + x, offset[1] + y, offset[2] + z);
1229 auto const density =
1230 m_boundary_density->get_node_value_at_boundary(global);
1231 Cell const local(x, y, z);
1232 ek::accessor::Scalar::set(density_field, density, local);
1233 }
1234 }) // WALBERLA_FOR_ALL_CELLS_XYZ
1235
1236 auto const bci = density_field->xyzSize();
1237 density_writer->set_content(
1239 density_writer->set_dims(
1240 Vector3<uint_t>(bci.xSize(), bci.ySize(), bci.zSize()));
1241 }
1242 };
1243
1244 vtk_obj.addBeforeFunction(std::move(before_function));
1245 vtk_obj.addCellDataWriter(density_writer);
1246 }
1247 if (flag_observables & static_cast<int>(EKOutputVTK::flux)) {
1248 auto const unit_conversion = FloatType_c(units.at("flux"));
1249 auto const blocks = m_lattice->get_blocks();
1253 auto before_function = [this, blocks, flux_writer]() {
1254 for (auto &block : *blocks) {
1256 auto const bci = flux_field->xyzSize();
1258 auto const block_offset = m_lattice->get_block_corner(block, true);
1259 std::size_t index = 0u;
1260 for (auto x = bci.xMin(); x <= bci.xMax(); ++x) {
1261 for (auto y = bci.yMin(); y <= bci.yMax(); ++y) {
1262 for (auto z = bci.zMin(); z <= bci.zMax(); ++z) {
1263 auto const node = block_offset + Utils::Vector3i{{x, y, z}};
1264 if (m_boundary_flux->node_is_boundary(node)) {
1265 auto const &vec =
1266 m_boundary_flux->get_node_value_at_boundary(node);
1267 values[3u * index + 0u] = vec[0];
1268 values[3u * index + 1u] = vec[1];
1269 values[3u * index + 2u] = vec[2];
1270 }
1271 ++index;
1272 }
1273 }
1274 }
1275 flux_writer->set_content(std::move(values));
1276 flux_writer->set_dims(
1277 Vector3<uint_t>(bci.xSize(), bci.ySize(), bci.zSize()));
1278 }
1279 };
1280 vtk_obj.addBeforeFunction(std::move(before_function));
1281 vtk_obj.addCellDataWriter(flux_writer);
1282 }
1283 if (flag_observables & static_cast<int>(EKOutputVTK::boundary)) {
1284 vtk_obj.addCellDataWriter(make_shared<BoundaryVTKWriter<float>>(
1286 }
1287 }
1288
1289 ~EKinWalberlaImpl() override = default;
1290};
1291
1292} // namespace walberla
Vector implementation and trait types for boost qvm interoperability.
Interface of a lattice-based electrokinetic 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
std::pair< Utils::Vector3i, Utils::Vector3i > get_local_grid_range(bool with_halo=false) const
Utils::Vector3i get_block_corner(IBlock const &block, bool lower) const
Boundary class optimized for sparse data.
vtk::BlockCellDataWriter< OutputType, 1u > Base
OutputType evaluate(cell_idx_t const x, cell_idx_t const y, cell_idx_t const z, cell_idx_t const) override
BoundaryVTKWriter(ConstBlockDataID const &flag_field_id, std::string const &id, FlagUID const &boundary_flag)
OutputType evaluate(cell_idx_t const x, cell_idx_t const y, cell_idx_t const z, cell_idx_t const) override
OutputType evaluate(cell_idx_t const x, cell_idx_t const y, cell_idx_t const z, cell_idx_t const f) override
void set_dims(Vector3< uint_t > dims)
std::size_t get_first_index(cell_idx_t const x, cell_idx_t const y, cell_idx_t const z)
VTKWriter(ConstBlockDataID const &block_id, std::string const &id, FloatType unit_conversion)
Class that runs and controls the EK on waLBerla.
~EKinWalberlaImpl() override=default
void set_slice_flux_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< std::optional< Utils::Vector3d > > const &flux) override
void set_kT(double kT) override
std::optional< double > get_node_density_at_boundary(Utils::Vector3i const &node, bool consider_ghosts=false) const override
double get_kT() const noexcept override
std::unique_ptr< DiffusiveFluxKernelElectrostatic > m_diffusive_flux_electrostatic
walberla::FlagField< walberla::uint8_t > FlagField
void set_friction_coupling(bool friction_coupling) override
void set_rng_state(uint64_t counter) override
void integrate(std::size_t potential_id, std::size_t velocity_id, std::size_t force_id, double lb_density) override
bool set_node_density_boundary(Utils::Vector3i const &node, double density) override
void set_slice_density(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< double > const &density) override
std::shared_ptr< FullCommunicator > m_full_communication
std::bitset< GhostComm::SIZE > m_pending_ghost_comm
std::vector< double > get_slice_density(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
void update_density_boundary_from_shape(const std::vector< int > &raster_flat, const std::vector< double > &data_flat) override
void update_flux_boundary_from_shape(const std::vector< int > &raster_flat, const std::vector< double > &data_flat) override
std::vector< bool > get_slice_is_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
void set_ext_efield(Utils::Vector3d const &field) override
bool set_node_density(Utils::Vector3i const &node, double density) override
double get_diffusion() const noexcept override
void set_valency(double valency) override
bool is_thermalized() const noexcept override
auto add_to_storage(std::string const tag, FloatType value)
Convenience function to add a field with a custom allocator.
std::optional< double > get_node_density(Utils::Vector3i const &node, bool consider_ghosts=false) const override
std::size_t get_density_id() const noexcept override
unsigned int get_seed() const noexcept override
std::unique_ptr< ContinuityKernel > m_continuity
double get_valency() const noexcept override
typename FieldTrait< FloatType, Architecture >::FluxField FluxField
std::shared_ptr< BoundaryModelFlux > m_boundary_flux
void set_slice_density_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< std::optional< double > > const &density) override
void reset_flux_boundary_handling(std::shared_ptr< BlockStorage > const &blocks)
bool get_advection() const noexcept override
std::optional< Utils::Vector3d > get_node_flux_vector(Utils::Vector3i const &node, bool consider_ghosts=false) const override
typename FieldTrait< FloatType, Architecture >::template BoundaryCommScheme< typename stencil::D3Q27 > BoundaryFullCommunicator
typename FieldTrait< FloatType, Architecture >::template RegularCommScheme< typename stencil::D3Q27 > FullCommunicator
std::optional< uint64_t > get_rng_state() const override
void set_advection(bool advection) override
bool remove_node_from_density_boundary(Utils::Vector3i const &node) override
typename FieldTrait< FloatType, Architecture >::template PackInfo< Field > PackInfo
bool get_friction_coupling() const noexcept override
Utils::Vector3d get_ext_efield() const noexcept override
stencil::D3Q27 Stencil
Stencil for collision and streaming operations.
LatticeWalberla const & get_lattice() const noexcept override
bool remove_node_from_flux_boundary(Utils::Vector3i const &node) override
bool is_double_precision() const noexcept override
typename FieldTrait< FloatType, Architecture >::DensityField DensityField
bool is_gpu() const noexcept override
void clear_density_boundaries() override
std::optional< bool > get_node_is_flux_boundary(Utils::Vector3i const &node, bool consider_ghosts) const override
LatticeWalberla::Lattice_T BlockStorage
Lattice model (e.g.
void reset_density_boundary_handling(std::shared_ptr< BlockStorage > const &blocks)
std::vector< std::optional< Utils::Vector3d > > get_slice_flux_at_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
ResourceObserver m_mpi_cart_comm_observer
std::vector< std::optional< double > > get_slice_density_at_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj, LatticeModel::units_map const &units, int flag_observables) override
std::optional< Utils::Vector3d > get_node_flux_at_boundary(Utils::Vector3i const &node, bool consider_ghosts=false) const override
std::vector< double > get_slice_flux_vector(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
typename FieldTrait< FloatType >::DensityField _DensityField
std::shared_ptr< BoundaryFullCommunicator > m_boundary_communicator
std::optional< bool > get_node_is_density_boundary(Utils::Vector3i const &node, bool consider_ghosts) const override
std::unique_ptr< DiffusiveFluxKernel > m_diffusive_flux
void set_diffusion(double diffusion) override
std::size_t stencil_size() const noexcept override
FlagUID const Boundary_flag
Flag for boundary cells.
bool set_node_flux_boundary(Utils::Vector3i const &node, Utils::Vector3d const &flux) override
typename FieldTrait< FloatType >::FluxField _FluxField
std::optional< bool > get_node_is_boundary(Utils::Vector3i const &node, bool consider_ghosts=false) const override
std::unique_ptr< BoundaryModelDensity > m_boundary_density
EKinWalberlaImpl(std::shared_ptr< LatticeWalberla > lattice, double diffusion, double kT, double valency, Utils::Vector3d const &ext_efield, double density, bool advection, bool friction_coupling, bool thermalized, unsigned int seed)
std::shared_ptr< LatticeWalberla > m_lattice
Block forest.
FlagUID const Domain_flag
Flag for domain cells, i.e.
void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override
void setup_boundary_handle(std::shared_ptr< LatticeWalberla > lattice, std::shared_ptr< Boundary_T > boundary)
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:175
STL namespace.
auto get_vector(GhostLayerField< double, uint_t{13u}> const *flux_field, Cell const &cell)
void initialize(GhostLayerField< double, uint_t{13u}> *flux_field, std::array< double, 13 > const &values)
void initialize(GhostLayerField< double, 1u > *scalar_field, double const &value)
void set(GhostLayerField< double, 1u > *scalar_field, double const &value, Cell const &cell)
auto get(GhostLayerField< double, 1u > const *scalar_field, Cell const &cell)
static FUNC_PREFIX double *RESTRICT const double *RESTRICT const int64_t const int64_t const int64_t const int64_t const int64_t const int64_t const int64_t const int64_t const int64_t const int64_t const uint32_t uint32_t uint32_t uint32_t uint32_t uint32_t uint32_t seed
\file PackInfoPdfDoublePrecision.cpp \author pystencils
auto to_vector3d(Vector3< T > const &v) noexcept
std::vector< double > fill_3D_scalar_array(std::vector< double > const &vec_flat, Utils::Vector3i const &grid_size)
Definition boundary.hpp:60
void set_boundary_from_grid(BoundaryModel &boundary, LatticeWalberla const &lattice, std::vector< int > const &raster_flat, std::vector< DataType > const &data_flat)
Definition boundary.hpp:81
void copy_block_buffer(CellInterval const &bci, CellInterval const &ci, Utils::Vector3i const &block_offset, Utils::Vector3i const &lower_corner, auto &&kernel)
Synchronize data between a sliced block and a container.
std::optional< BlockAndCell > get_block_and_cell(::LatticeWalberla const &lattice, signed_integral_vector auto const &node, bool consider_ghost_layers)
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.
std::optional< walberla::cell::CellInterval > get_block_interval(::LatticeWalberla const &lattice, Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, Utils::Vector3i const &block_offset, IBlock const &block)
std::optional< walberla::cell::CellInterval > get_interval(::LatticeWalberla const &lattice, Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner)
std::vector< Utils::Vector3d > fill_3D_vector_array(std::vector< double > const &vec_flat, Utils::Vector3i const &grid_size)
Definition boundary.hpp:36
Observer to monitor the lifetime of a shared resource.
blockforest::communication::UniformBufferedScheme< Stencil > RegularCommScheme
GhostLayerField< FT, FluxCount > FluxField
field::communication::PackInfo< Field > PackInfo
blockforest::communication::UniformBufferedScheme< Stencil > BoundaryCommScheme
GhostCommFlags
Ghost communication operations.