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/vtk/FlagFieldCellFilter.h>
29#include <field/vtk/VTKWriter.h>
30#include <stencil/D3Q27.h>
31#include <waLBerlaDefinitions.h>
32#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
33#include <gpu/AddGPUFieldToStorage.h>
34#include <gpu/communication/MemcpyPackInfo.h>
35#include <gpu/communication/UniformGPUScheme.h>
36#endif
37
38#include "../BoundaryHandling.hpp"
39#include "../BoundaryPackInfo.hpp"
40#include "../utils/boundary.hpp"
41#include "../utils/types_conversion.hpp"
42#include "ek_kernels.hpp"
43#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
44#include "ek_kernels.cuh"
45#endif
46
53
54#include <utils/Vector.hpp>
55
56#include <cstddef>
57#include <iterator>
58#include <memory>
59#include <optional>
60#include <stdexcept>
61#include <string>
62#include <type_traits>
63#include <variant>
64#include <vector>
65
66namespace walberla {
67
68/** @brief Class that runs and controls the EK on waLBerla. */
69template <std::size_t FluxCount = 13, typename FloatType = double,
70 lbmpy::Arch Architecture = lbmpy::Arch::CPU>
72#if not defined(WALBERLA_BUILD_WITH_CUDA)
73 static_assert(Architecture != lbmpy::Arch::GPU,
74 "waLBerla was compiled without CUDA support");
75#endif
76 using ContinuityKernel =
78 using DiffusiveFluxKernelUnthermalized =
79 typename detail::KernelTrait<FloatType,
80 Architecture>::DiffusiveFluxKernel;
81 using DiffusiveFluxKernelThermalized = typename detail::KernelTrait<
82 FloatType, Architecture>::DiffusiveFluxKernelThermalized;
83 using AdvectiveFluxKernel =
84 typename detail::KernelTrait<FloatType,
85 Architecture>::AdvectiveFluxKernel;
86 using FrictionCouplingKernel =
87 typename detail::KernelTrait<FloatType,
88 Architecture>::FrictionCouplingKernel;
89 using DiffusiveFluxKernelElectrostaticUnthermalized =
90 typename detail::KernelTrait<
91 FloatType, Architecture>::DiffusiveFluxKernelElectrostatic;
92 using DiffusiveFluxKernelElectrostaticThermalized =
93 typename detail::KernelTrait<
94 FloatType, Architecture>::DiffusiveFluxKernelElectrostaticThermalized;
95
96 using DiffusiveFluxKernel = std::variant<DiffusiveFluxKernelUnthermalized,
97 DiffusiveFluxKernelThermalized>;
98 using DiffusiveFluxKernelElectrostatic =
99 std::variant<DiffusiveFluxKernelElectrostaticUnthermalized,
100 DiffusiveFluxKernelElectrostaticThermalized>;
101
102 using Dirichlet =
104 using FixedFlux =
106
109 using BoundaryModelFlux =
111
112public:
113 /** @brief Stencil for collision and streaming operations. */
114 using Stencil = stencil::D3Q27;
115 /** @brief Lattice model (e.g. blockforest). */
117
118protected:
119 template <typename FT, lbmpy::Arch AT = lbmpy::Arch::CPU> struct FieldTrait {
120 // Type definitions
121 using FluxField = GhostLayerField<FT, FluxCount>;
122 using DensityField = GhostLayerField<FT, 1>;
123 template <class Field>
124 using PackInfo = field::communication::PackInfo<Field>;
125 template <class Stencil>
127 blockforest::communication::UniformBufferedScheme<Stencil>;
128 template <class Stencil>
130 blockforest::communication::UniformBufferedScheme<Stencil>;
131 };
132 using FlagField = walberla::FlagField<walberla::uint8_t>;
133#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
134 template <typename FT> struct FieldTrait<FT, lbmpy::Arch::GPU> {
135 private:
136 static auto constexpr AT = lbmpy::Arch::GPU;
137 template <class Field>
138 using MemcpyPackInfo = gpu::communication::MemcpyPackInfo<Field>;
139
140 public:
141 template <typename Stencil>
142 class UniformGPUScheme
143 : public gpu::communication::UniformGPUScheme<Stencil> {
144 public:
145 explicit UniformGPUScheme(auto const &bf)
146 : gpu::communication::UniformGPUScheme<Stencil>(
147 bf, /* sendDirectlyFromGPU */ false,
148 /* useLocalCommunication */ false) {}
149 };
150 using FluxField = gpu::GPUField<FT>;
151 using DensityField = gpu::GPUField<FT>;
152 template <class Field> using PackInfo = MemcpyPackInfo<Field>;
153 template <class Stencil>
154 using RegularCommScheme = UniformGPUScheme<Stencil>;
155 template <class Stencil>
156 using BoundaryCommScheme =
157 blockforest::communication::UniformBufferedScheme<Stencil>;
158 };
159 using GPUField = gpu::GPUField<FloatType>;
160#endif
161
162 struct GhostComm {
163 /** @brief Ghost communication operations. */
164 enum GhostCommFlags : unsigned {
165 FLB, ///< flux boundary communication
166 DENS, ///< density communication
167 SIZE
168 };
169 };
170
171 // "underlying" field types (`GPUField` has no f-size info at compile time)
174
175public:
179
180 template <typename T> FloatType FloatType_c(T t) {
181 return numeric_cast<FloatType>(t);
182 }
183
184 [[nodiscard]] std::size_t stencil_size() const noexcept override {
185 return FluxCount;
186 }
187
188 [[nodiscard]] bool is_double_precision() const noexcept override {
189 return std::is_same_v<FloatType, double>;
190 }
191
192private:
193 FloatType m_diffusion;
194 FloatType m_kT;
195 FloatType m_valency;
196 Utils::Vector3d m_ext_efield;
197 bool m_advection;
198 bool m_friction_coupling;
199 unsigned int m_seed;
200
201protected:
202 // Block data access handles
204
205 BlockDataID m_flux_field_id;
206
209
210 /** Flag for domain cells, i.e. all cells. */
211 FlagUID const Domain_flag{"domain"};
212 /** Flag for boundary cells. */
213 FlagUID const Boundary_flag{"boundary"};
214
215 /** Block forest */
216 std::shared_ptr<LatticeWalberla> m_lattice;
217
218 std::unique_ptr<BoundaryModelDensity> m_boundary_density;
219 std::shared_ptr<BoundaryModelFlux> m_boundary_flux;
220
221 std::unique_ptr<DiffusiveFluxKernel> m_diffusive_flux;
222 std::unique_ptr<DiffusiveFluxKernelElectrostatic>
224 std::unique_ptr<ContinuityKernel> m_continuity;
225
226 // ResetFlux + external force
227 // TODO: kernel for that
228 // std::shared_ptr<ResetForce<PdfField, VectorField>> m_reset_force;
229
230 /**
231 * @brief Convenience function to add a field with a custom allocator.
232 *
233 * When vectorization is off, let waLBerla decide which memory allocator
234 * to use. When vectorization is on, the aligned memory allocator is
235 * required, otherwise <tt>cpu_vectorize_info["assume_aligned"]</tt> will
236 * trigger assertions. That is because for single-precision kernels the
237 * waLBerla heuristic in <tt>src/field/allocation/FieldAllocator.h</tt>
238 * will fall back to @c StdFieldAlloc, yet @c AllocateAligned is needed
239 * for intrinsics to work.
240 */
241 template <typename Field>
242 auto add_to_storage(std::string const tag, FloatType value) {
243 auto const &blocks = m_lattice->get_blocks();
244 auto const n_ghost_layers = m_lattice->get_ghost_layers();
245 if constexpr (Architecture == lbmpy::Arch::CPU) {
246 return field::addToStorage<Field>(blocks, tag, FloatType{value},
247 field::fzyx, n_ghost_layers);
248 }
249#if defined(__CUDACC__) and defined(WALBERLA_BUILD_WITH_CUDA)
250 else {
251 auto field_id = gpu::addGPUFieldToStorage<GPUField>(
252 blocks, tag, Field::F_SIZE, field::fzyx, n_ghost_layers);
253 if constexpr (std::is_same_v<Field, _DensityField>) {
254 for (auto block = blocks->begin(); block != blocks->end(); ++block) {
255 auto field = block->template getData<GPUField>(field_id);
256 ek::accessor::Scalar::initialize(field, FloatType{value});
257 }
258 } else if constexpr (std::is_same_v<Field, _FluxField>) {
259 for (auto block = blocks->begin(); block != blocks->end(); ++block) {
260 auto field = block->template getData<GPUField>(field_id);
262 std::array<FloatType, FluxCount>{});
263 }
264 }
265 return field_id;
266 }
267#endif
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
315 add_to_storage<_DensityField>("density field", FloatType_c(density));
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<
348 flux_boundary_packinfo->setup_boundary_handle(m_lattice, m_boundary_flux);
349 m_boundary_communicator->addPackInfo(flux_boundary_packinfo);
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 }
365 [[nodiscard]] bool get_friction_coupling() const noexcept override {
366 return m_friction_coupling;
367 }
368 [[nodiscard]] Utils::Vector3d get_ext_efield() const noexcept override {
369 return m_ext_efield;
370 }
371 [[nodiscard]] bool is_thermalized() const noexcept override {
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);
393 std::visit(visitor, *m_diffusive_flux_electrostatic);
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
411 void set_friction_coupling(bool friction_coupling) override {
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(
467 m_flux_field_id, BlockDataID{}, m_density_field_id,
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(
486 m_flux_field_id, BlockDataID{}, m_density_field_id,
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
533 auto *kernel_electrostatic =
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,
543 BlockDataID(velocity_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(
552 BlockDataID(force_id), m_flux_field_id, FloatType_c(get_diffusion()),
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 &it : m_vtk_auto) {
591 auto &vtk_handle = it.second;
592 if (vtk_handle->enabled) {
593 vtk::writeFiles(vtk_handle->ptr)();
594 vtk_handle->execution_count++;
595 }
596 }
597 }
598
599public:
600 void integrate(std::size_t potential_id, std::size_t velocity_id,
601 std::size_t force_id, double lb_density) override {
602
604 update_boundary_fields();
605
606 if (get_diffusion() == 0.)
607 return;
608
609 if (get_valency() != 0.) {
610 if (potential_id == walberla::BlockDataID{}) {
611 throw std::runtime_error("Walberla EK: electrostatic potential enabled "
612 "but no field accessible. potential id is " +
613 std::to_string(potential_id));
614 }
615 kernel_diffusion_electrostatic(potential_id);
616 } else {
617 kernel_diffusion();
618 }
619
620 kernel_migration();
621 kernel_boundary_flux();
622 // friction coupling
623 if (get_friction_coupling()) {
624 if (force_id == walberla::BlockDataID{}) {
625 throw std::runtime_error("Walberla EK: friction coupling enabled but "
626 "no force field accessible. force_id is " +
627 std::to_string(force_id) +
628 ". Hint: LB may be inactive.");
629 }
630 kernel_friction_coupling(force_id, lb_density);
631 }
632
633 if (get_advection()) {
634 if (velocity_id == walberla::BlockDataID{}) {
635 throw std::runtime_error("Walberla EK: advection enabled but no "
636 "velocity field accessible. velocity_id is " +
637 std::to_string(velocity_id) +
638 ". Hint: LB may be inactive.");
639 }
640 kernel_advection(velocity_id);
641 kernel_boundary_flux();
642 }
643 kernel_continuity();
644
645 // is this the expected behavior when reactions are included?
646 kernel_boundary_density();
648
649 // Handle VTK writers
651 }
652
653 [[nodiscard]] std::size_t get_density_id() const noexcept override {
654 static_assert(std::is_same_v<std::size_t, walberla::uint_t>);
655 return static_cast<std::size_t>(m_density_field_id);
656 }
657
658 bool set_node_density(Utils::Vector3i const &node, double density) override {
660 auto bc = get_block_and_cell(get_lattice(), node, false);
661 if (!bc)
662 return false;
663
664 auto density_field =
665 bc->block->template getData<DensityField>(m_density_field_id);
666 ek::accessor::Scalar::set(density_field, FloatType_c(density), bc->cell);
667 return true;
668 }
669
670 [[nodiscard]] std::optional<double>
672 bool consider_ghosts = false) const override {
673 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
674
675 if (!bc)
676 return std::nullopt;
677
678 auto const density_field =
679 bc->block->template getData<DensityField>(m_density_field_id);
680 return {double_c(ek::accessor::Scalar::get(density_field, bc->cell))};
681 }
682
683 [[nodiscard]] std::vector<double>
685 Utils::Vector3i const &upper_corner) const override {
686 std::vector<double> out;
687#ifndef NDEBUG
688 uint_t values_size = 0u;
689#endif
690 auto const &lattice = get_lattice();
691 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
692 out = std::vector<double>(ci->numCells());
693 for (auto &block : *lattice.get_blocks()) {
694 auto const block_offset = lattice.get_block_corner(block, true);
695 if (auto const bci = get_block_interval(
696 lattice, lower_corner, upper_corner, block_offset, block)) {
697 auto const density_field =
698 block.template getData<DensityField>(m_density_field_id);
699 auto const values = ek::accessor::Scalar::get(density_field, *bci);
700 assert(values.size() == bci->numCells());
701#ifndef NDEBUG
702 values_size += bci->numCells();
703#endif
704 auto kernel = [&values, &out](unsigned const block_index,
705 unsigned const local_index,
706 Utils::Vector3i const &) {
707 out[local_index] = double_c(values[block_index]);
708 };
709
710 copy_block_buffer(*bci, *ci, block_offset, lower_corner, kernel);
711 }
712 }
713 assert(values_size == ci->numCells());
714 }
715 return out;
716 }
717
718 void set_slice_density(Utils::Vector3i const &lower_corner,
719 Utils::Vector3i const &upper_corner,
720 std::vector<double> const &density) override {
722 auto const &lattice = get_lattice();
723 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
724 assert(density.size() == ci->numCells());
725 for (auto &block : *lattice.get_blocks()) {
726 auto const block_offset = lattice.get_block_corner(block, true);
727 if (auto const bci = get_block_interval(
728 lattice, lower_corner, upper_corner, block_offset, block)) {
729 auto const density_field =
730 block.template getData<DensityField>(m_density_field_id);
731 std::vector<FloatType> values(bci->numCells());
732
733 auto kernel = [&values, &density](unsigned const block_index,
734 unsigned const local_index,
735 Utils::Vector3i const &) {
736 values[block_index] = numeric_cast<FloatType>(density[local_index]);
737 };
738
739 copy_block_buffer(*bci, *ci, block_offset, lower_corner, kernel);
740 ek::accessor::Scalar::set(density_field, values, *bci);
741 }
742 }
743 }
744 }
745
746 [[nodiscard]] std::optional<Utils::Vector3d>
748 bool consider_ghosts = false) const override {
749 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
750
751 if (!bc)
752 return std::nullopt;
753
754 auto const flux_field =
755 bc->block->template getData<FluxField>(m_flux_field_id);
756 return to_vector3d(ek::accessor::Flux::get_vector(flux_field, bc->cell));
757 }
758
759 std::vector<double>
761 Utils::Vector3i const &upper_corner) const override {
762 std::vector<double> out;
763#ifndef NDEBUG
764 uint_t values_size = 0;
765#endif
766 auto const &lattice = get_lattice();
767 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
768 out = std::vector<double>(3u * ci->numCells());
769 for (auto &block : *lattice.get_blocks()) {
770 auto const block_offset = lattice.get_block_corner(block, true);
771 if (auto const bci = get_block_interval(
772 lattice, lower_corner, upper_corner, block_offset, block)) {
773 auto const flux_field =
774 block.template getData<FluxField>(m_flux_field_id);
775 auto const values = ek::accessor::Flux::get_vector(flux_field, *bci);
776 assert(values.size() == 3u * bci->numCells());
777#ifndef NDEBUG
778 values_size += 3u * bci->numCells();
779#endif
780
781 auto kernel = [&values, &out, this](unsigned const block_index,
782 unsigned const local_index,
783 Utils::Vector3i const &node) {
784 if (m_boundary_flux->node_is_boundary(node)) {
785 auto const &vec =
786 m_boundary_flux->get_node_value_at_boundary(node);
787 for (uint_t f = 0u; f < 3u; ++f) {
788 out[3u * local_index + f] = double_c(vec[f]);
789 }
790 } else {
791 for (uint_t f = 0u; f < 3u; ++f) {
792 out[3u * local_index + f] =
793 double_c(values[3u * block_index + f]);
794 }
795 }
796 };
797
798 copy_block_buffer(*bci, *ci, block_offset, lower_corner, kernel);
799 }
800 }
801 assert(values_size == 3u * ci->numCells());
802 }
803 return out;
804 }
805
810
811 void clear_density_boundaries() override {
813 }
814
816 Utils::Vector3d const &flux) override {
818 auto bc = get_block_and_cell(get_lattice(), node, true);
819 if (!bc)
820 return false;
821
822 m_boundary_flux->set_node_value_at_boundary(
823 node, to_vector3<FloatType>(flux), *bc);
824 return true;
825 }
826
827 [[nodiscard]] std::optional<Utils::Vector3d>
829 bool consider_ghosts = false) const override {
830 assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::FLB)));
831 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
832 if (!bc or !m_boundary_flux->node_is_boundary(node))
833 return std::nullopt;
834
835 return {to_vector3d(m_boundary_flux->get_node_value_at_boundary(node))};
836 }
837
840 auto bc = get_block_and_cell(get_lattice(), node, true);
841 if (!bc)
842 return false;
843
844 m_boundary_flux->remove_node_from_boundary(node, *bc);
845 return true;
846 }
847
849 double density) override {
850 auto bc = get_block_and_cell(get_lattice(), node, true);
851 if (!bc)
852 return false;
853
854 m_boundary_density->set_node_value_at_boundary(node, FloatType_c(density),
855 *bc);
856
857 return true;
858 }
859
860 [[nodiscard]] std::optional<double>
862 bool consider_ghosts = false) const override {
863 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
864 if (!bc or !m_boundary_density->node_is_boundary(node))
865 return std::nullopt;
866
867 return {double_c(m_boundary_density->get_node_value_at_boundary(node))};
868 }
869
871 Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
872 std::vector<std::optional<double>> const &density) override {
873 auto const &lattice = get_lattice();
874 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
875 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
876 auto const lower_cell = ci->min();
877 auto const upper_cell = ci->max();
878 auto it = density.begin();
879 assert(density.size() == ci->numCells());
880 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
881 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
882 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
883 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
884 auto const bc = get_block_and_cell(lattice, node, false);
885 auto const &opt = *it;
886 if (opt) {
887 m_boundary_density->set_node_value_at_boundary(
888 node, FloatType_c(*opt), *bc);
889 } else {
890 m_boundary_density->remove_node_from_boundary(node, *bc);
891 }
892 ++it;
893 }
894 }
895 }
896 }
897 }
898
899 [[nodiscard]] std::vector<std::optional<double>>
901 Utils::Vector3i const &lower_corner,
902 Utils::Vector3i const &upper_corner) const override {
903 std::vector<std::optional<double>> out;
904 auto const &lattice = get_lattice();
905 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
906 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
907 auto const lower_cell = ci->min();
908 auto const upper_cell = ci->max();
909 auto const n_values = ci->numCells();
910 out.reserve(n_values);
911 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
912 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
913 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
914 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
915 if (m_boundary_density->node_is_boundary(node)) {
916 out.emplace_back(double_c(
917 m_boundary_density->get_node_value_at_boundary(node)));
918 } else {
919 out.emplace_back(std::nullopt);
920 }
921 }
922 }
923 }
924 assert(out.size() == n_values);
925 }
926 return out;
927 }
928
930 Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
931 std::vector<std::optional<Utils::Vector3d>> const &flux) override {
933 auto const &lattice = get_lattice();
934 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
935 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
936 auto const lower_cell = ci->min();
937 auto const upper_cell = ci->max();
938 auto it = flux.begin();
939 assert(flux.size() == ci->numCells());
940 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
941 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
942 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
943 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
944 auto const bc = get_block_and_cell(lattice, node, false);
945 auto const &opt = *it;
946 if (opt) {
947 m_boundary_flux->set_node_value_at_boundary(
948 node, to_vector3<FloatType>(*opt), *bc);
949 } else {
950 m_boundary_flux->remove_node_from_boundary(node, *bc);
951 }
952 ++it;
953 }
954 }
955 }
956 }
957 }
958
959 [[nodiscard]] std::vector<std::optional<Utils::Vector3d>>
961 Utils::Vector3i const &lower_corner,
962 Utils::Vector3i const &upper_corner) const override {
963 std::vector<std::optional<Utils::Vector3d>> out;
964 auto const &lattice = get_lattice();
965 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
966 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
967 auto const lower_cell = ci->min();
968 auto const upper_cell = ci->max();
969 auto const n_values = ci->numCells();
970 out.reserve(n_values);
971 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
972 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
973 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
974 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
975 if (m_boundary_flux->node_is_boundary(node)) {
976 out.emplace_back(to_vector3d(
977 m_boundary_flux->get_node_value_at_boundary(node)));
978 } else {
979 out.emplace_back(std::nullopt);
980 }
981 }
982 }
983 }
984 assert(out.size() == n_values);
985 }
986 return out;
987 }
988
989 [[nodiscard]] std::vector<bool>
991 Utils::Vector3i const &upper_corner) const override {
992 std::vector<bool> out;
993 auto const &lattice = get_lattice();
994 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
995 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
996 auto const lower_cell = ci->min();
997 auto const upper_cell = ci->max();
998 auto const n_values = ci->numCells();
999 out.reserve(n_values);
1000 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
1001 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
1002 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
1003 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
1004 out.emplace_back(m_boundary_density->node_is_boundary(node) or
1005 m_boundary_flux->node_is_boundary(node));
1006 }
1007 }
1008 }
1009 assert(out.size() == n_values);
1010 }
1011 return out;
1012 }
1013
1015 auto bc = get_block_and_cell(get_lattice(), node, true);
1016 if (!bc)
1017 return false;
1018
1019 m_boundary_density->remove_node_from_boundary(node, *bc);
1020
1021 return true;
1022 }
1023
1024 [[nodiscard]] std::optional<bool>
1026 bool consider_ghosts) const override {
1027 assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::FLB)));
1028 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
1029 if (!bc)
1030 return std::nullopt;
1031
1032 return {m_boundary_flux->node_is_boundary(node)};
1033 }
1034
1035 [[nodiscard]] std::optional<bool>
1037 bool consider_ghosts) const override {
1038 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
1039 if (!bc)
1040 return std::nullopt;
1041
1042 return {m_boundary_density->node_is_boundary(node)};
1043 }
1044
1045 [[nodiscard]] std::optional<bool>
1047 bool consider_ghosts = false) const override {
1048 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
1049 if (!bc)
1050 return std::nullopt;
1051
1052 return {m_boundary_density->node_is_boundary(node) or
1053 m_boundary_flux->node_is_boundary(node)};
1054 }
1055
1057 const std::vector<int> &raster_flat,
1058 const std::vector<double> &data_flat) override {
1060 auto const grid_size = get_lattice().get_grid_dimensions();
1061 auto const data = fill_3D_vector_array(data_flat, grid_size);
1062 set_boundary_from_grid(*m_boundary_flux, get_lattice(), raster_flat, data);
1064 }
1065
1067 const std::vector<int> &raster_flat,
1068 const std::vector<double> &data_flat) override {
1069 auto const grid_size = get_lattice().get_grid_dimensions();
1070 auto const data = fill_3D_scalar_array(data_flat, grid_size);
1072 data);
1074 }
1075
1077
1079 m_boundary_density->boundary_update();
1080 }
1081
1082 [[nodiscard]] LatticeWalberla const &get_lattice() const noexcept override {
1083 return *m_lattice;
1084 }
1085
1086 [[nodiscard]] bool is_gpu() const noexcept override {
1087 return Architecture == lbmpy::Arch::GPU;
1088 }
1089
1090 void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override {
1091 field::FlagFieldCellFilter<FlagField> fluid_filter(m_flag_field_density_id);
1092 fluid_filter.addFlag(Boundary_flag);
1093 vtk_obj.addCellExclusionFilter(fluid_filter);
1094 }
1095
1096protected:
1097 template <typename VecType, uint_t F_SIZE_ARG, typename OutputType>
1098 class VTKWriter : public vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG> {
1099 public:
1100 VTKWriter(ConstBlockDataID const &block_id, std::string const &id,
1101 FloatType unit_conversion)
1102 : vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG>(id),
1103 m_conversion(unit_conversion), m_content{} {}
1104
1105 protected:
1106 void configure() override { WALBERLA_ASSERT_NOT_NULLPTR(this->block_); }
1107
1108 std::size_t get_first_index(cell_idx_t const x, cell_idx_t const y,
1109 cell_idx_t const z) {
1110 return (static_cast<std::size_t>(x) * m_dims[2] * m_dims[1] +
1111 static_cast<std::size_t>(y) * m_dims[2] +
1112 static_cast<std::size_t>(z)) *
1113 F_SIZE_ARG;
1114 }
1115
1116 FloatType m_conversion;
1117 VecType m_content;
1118 Vector3<uint_t> m_dims;
1119
1120 public:
1121 void set_content(VecType content) { m_content = std::move(content); }
1122
1123 void set_dims(Vector3<uint_t> dims) { m_dims = dims; }
1124 };
1125
1126 template <typename OutputType = float>
1128 : public VTKWriter<std::vector<FloatType>, 1u, OutputType> {
1129 public:
1130 using Base = VTKWriter<std::vector<FloatType>, 1u, OutputType>;
1131 using Base::Base;
1132 using Base::evaluate;
1133
1134 protected:
1135 OutputType evaluate(cell_idx_t const x, cell_idx_t const y,
1136 cell_idx_t const z, cell_idx_t const) override {
1137 WALBERLA_ASSERT(!this->m_content.empty());
1138 auto const density = this->m_content[this->get_first_index(x, y, z)];
1139 return numeric_cast<OutputType>(this->m_conversion * density);
1140 }
1141 };
1142
1143 template <typename OutputType = float>
1145 : public VTKWriter<std::vector<FloatType>, 3u, OutputType> {
1146 public:
1147 using Base = VTKWriter<std::vector<FloatType>, 3u, OutputType>;
1148 using Base::Base;
1149 using Base::evaluate;
1150
1151 protected:
1152 OutputType evaluate(cell_idx_t const x, cell_idx_t const y,
1153 cell_idx_t const z, cell_idx_t const f) override {
1154 WALBERLA_ASSERT(!this->m_content.empty());
1155 auto flux = this->m_content[this->get_first_index(x, y, z) + f];
1156 return numeric_cast<OutputType>(this->m_conversion * flux);
1157 }
1158 };
1159
1160public:
1161 void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj,
1162 LatticeModel::units_map const &units,
1163 int flag_observables) override {
1164 if (flag_observables & static_cast<int>(EKOutputVTK::density)) {
1165 auto const unit_conversion = FloatType_c(units.at("density"));
1166 auto const blocks = m_lattice->get_blocks();
1167 WALBERLA_ASSERT_NOT_NULLPTR(blocks);
1168 auto density_writer = make_shared<DensityVTKWriter<float>>(
1169 m_density_field_id, "density", unit_conversion);
1170 auto before_function = [this, blocks, density_writer]() {
1171 for (auto &block : *blocks) {
1172 auto *density_field =
1173 block.template getData<DensityField>(m_density_field_id);
1174 auto const bci = density_field->xyzSize();
1175 density_writer->set_content(
1176 ek::accessor::Scalar::get(density_field, bci));
1177 density_writer->set_dims(Vector3<uint_t>(
1178 uint_c(bci.xSize()), uint_c(bci.ySize()), uint_c(bci.zSize())));
1179 }
1180 };
1181 vtk_obj.addBeforeFunction(std::move(before_function));
1182 vtk_obj.addCellDataWriter(density_writer);
1183 }
1184 if (flag_observables & static_cast<int>(EKOutputVTK::flux)) {
1185 auto const unit_conversion = FloatType_c(units.at("flux"));
1186 auto const blocks = m_lattice->get_blocks();
1187 WALBERLA_ASSERT_NOT_NULLPTR(blocks);
1188 auto flux_writer = make_shared<FluxVTKWriter<float>>(
1189 m_flux_field_id, "flux", unit_conversion);
1190 auto before_function = [this, blocks, flux_writer]() {
1191 for (auto &block : *blocks) {
1192 auto *flux_field = block.template getData<FluxField>(m_flux_field_id);
1193 auto const bci = flux_field->xyzSize();
1194 auto values = ek::accessor::Flux::get_vector(flux_field, bci);
1195 auto const block_offset = m_lattice->get_block_corner(block, true);
1196 std::size_t index = 0u;
1197 for (auto x = bci.xMin(); x <= bci.xMax(); ++x) {
1198 for (auto y = bci.yMin(); y <= bci.yMax(); ++y) {
1199 for (auto z = bci.zMin(); z <= bci.zMax(); ++z) {
1200 auto const node = block_offset + Utils::Vector3i{{x, y, z}};
1201 if (m_boundary_flux->node_is_boundary(node)) {
1202 auto const &vec =
1203 m_boundary_flux->get_node_value_at_boundary(node);
1204 values[3u * index + 0u] = vec[0];
1205 values[3u * index + 1u] = vec[1];
1206 values[3u * index + 2u] = vec[2];
1207 }
1208 ++index;
1209 }
1210 }
1211 }
1212 flux_writer->set_content(std::move(values));
1213 flux_writer->set_dims(Vector3<uint_t>(
1214 uint_c(bci.xSize()), uint_c(bci.ySize()), uint_c(bci.zSize())));
1215 }
1216 };
1217 vtk_obj.addBeforeFunction(std::move(before_function));
1218 vtk_obj.addCellDataWriter(flux_writer);
1219 }
1220 }
1221
1222 ~EKinWalberlaImpl() override = default;
1223};
1224
1225} // 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.
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)
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.