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-2023 The ESPResSo project
3 *
4 * This file is part of ESPResSo.
5 *
6 * ESPResSo is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation, either version 3 of the License, or
9 * (at your option) any later version.
10 *
11 * ESPResSo is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with this program. If not, see <http://www.gnu.org/licenses/>.
18 */
19
20#pragma once
21
22#include <blockforest/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#if defined(__CUDACC__)
32#include <gpu/AddGPUFieldToStorage.h>
33#include <gpu/HostFieldAllocator.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__)
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 using ContinuityKernel =
74 using DiffusiveFluxKernelUnthermalized =
75 typename detail::KernelTrait<FloatType,
76 Architecture>::DiffusiveFluxKernel;
77 using DiffusiveFluxKernelThermalized = typename detail::KernelTrait<
78 FloatType, Architecture>::DiffusiveFluxKernelThermalized;
79 using AdvectiveFluxKernel =
80 typename detail::KernelTrait<FloatType,
81 Architecture>::AdvectiveFluxKernel;
82 using FrictionCouplingKernel =
83 typename detail::KernelTrait<FloatType,
84 Architecture>::FrictionCouplingKernel;
85 using DiffusiveFluxKernelElectrostaticUnthermalized =
86 typename detail::KernelTrait<
87 FloatType, Architecture>::DiffusiveFluxKernelElectrostatic;
88 using DiffusiveFluxKernelElectrostaticThermalized =
89 typename detail::KernelTrait<
90 FloatType, Architecture>::DiffusiveFluxKernelElectrostaticThermalized;
91
92 using DiffusiveFluxKernel = std::variant<DiffusiveFluxKernelUnthermalized,
93 DiffusiveFluxKernelThermalized>;
94 using DiffusiveFluxKernelElectrostatic =
95 std::variant<DiffusiveFluxKernelElectrostaticUnthermalized,
96 DiffusiveFluxKernelElectrostaticThermalized>;
97
98 using Dirichlet =
100 using FixedFlux =
102
105 using BoundaryModelFlux =
107
108public:
109 /** @brief Stencil for collision and streaming operations. */
110 using Stencil = stencil::D3Q27;
111 /** @brief Lattice model (e.g. blockforest). */
113
114protected:
115 template <typename FT, lbmpy::Arch AT = lbmpy::Arch::CPU> struct FieldTrait {
116 // Type definitions
117 using FluxField = GhostLayerField<FT, FluxCount>;
118 using DensityField = GhostLayerField<FT, 1>;
119 template <class Field>
120 using PackInfo = field::communication::PackInfo<Field>;
121 template <class Stencil>
123 blockforest::communication::UniformBufferedScheme<Stencil>;
124 template <class Stencil>
126 blockforest::communication::UniformBufferedScheme<Stencil>;
127 };
128 using FlagField = walberla::FlagField<walberla::uint8_t>;
129#if defined(__CUDACC__)
130 template <typename FT> struct FieldTrait<FT, lbmpy::Arch::GPU> {
131 private:
132 static auto constexpr AT = lbmpy::Arch::GPU;
133 template <class Field>
134 using MemcpyPackInfo = gpu::communication::MemcpyPackInfo<Field>;
135
136 public:
137 template <typename Stencil>
138 class UniformGPUScheme
139 : public gpu::communication::UniformGPUScheme<Stencil> {
140 public:
141 explicit UniformGPUScheme(auto const &bf)
142 : gpu::communication::UniformGPUScheme<Stencil>(
143 bf, /* sendDirectlyFromGPU */ false,
144 /* useLocalCommunication */ false) {}
145 };
146 using FluxField = gpu::GPUField<FT>;
147 using DensityField = gpu::GPUField<FT>;
148 template <class Field> using PackInfo = MemcpyPackInfo<Field>;
149 template <class Stencil>
150 using RegularCommScheme = UniformGPUScheme<Stencil>;
151 template <class Stencil>
152 using BoundaryCommScheme =
153 blockforest::communication::UniformBufferedScheme<Stencil>;
154 };
155 using GPUField = gpu::GPUField<FloatType>;
156#endif
157
158 struct GhostComm {
159 /** @brief Ghost communication operations. */
160 enum GhostCommFlags : unsigned {
161 FLB, ///< flux boundary communication
162 DENS, ///< density communication
163 SIZE
164 };
165 };
166
167 // "underlying" field types (`GPUField` has no f-size info at compile time)
170
171public:
175
176#if defined(__CUDACC__)
179#endif
180
181 template <typename T> FloatType FloatType_c(T t) {
182 return numeric_cast<FloatType>(t);
183 }
184
185 [[nodiscard]] std::size_t stencil_size() const noexcept override {
186 return FluxCount;
187 }
188
189 [[nodiscard]] bool is_double_precision() const noexcept override {
190 return std::is_same_v<FloatType, double>;
191 }
192
193private:
194 FloatType m_diffusion;
195 FloatType m_kT;
196 FloatType m_valency;
197 Utils::Vector3d m_ext_efield;
198 bool m_advection;
199 bool m_friction_coupling;
200 unsigned int m_seed;
201
202protected:
203 // Block data access handles
205
206 BlockDataID m_flux_field_id;
207
210
211#if defined(__CUDACC__)
212 std::optional<BlockDataID> m_density_cpu_field_id;
213 std::optional<BlockDataID> m_flux_cpu_field_id;
214#endif
215
216 /** Flag for domain cells, i.e. all cells. */
217 FlagUID const Domain_flag{"domain"};
218 /** Flag for boundary cells. */
219 FlagUID const Boundary_flag{"boundary"};
220
221 /** Block forest */
222 std::shared_ptr<LatticeWalberla> m_lattice;
223
224 std::unique_ptr<BoundaryModelDensity> m_boundary_density;
225 std::shared_ptr<BoundaryModelFlux> m_boundary_flux;
226
227 std::unique_ptr<DiffusiveFluxKernel> m_diffusive_flux;
228 std::unique_ptr<DiffusiveFluxKernelElectrostatic>
230 std::unique_ptr<ContinuityKernel> m_continuity;
231
232#if defined(__CUDACC__)
233 std::shared_ptr<gpu::HostFieldAllocator<FloatType>> m_host_field_allocator;
234#endif
235
236 // ResetFlux + external force
237 // TODO: kernel for that
238 // std::shared_ptr<ResetForce<PdfField, VectorField>> m_reset_force;
239
240 /**
241 * @brief Convenience function to add a field with a custom allocator.
242 *
243 * When vectorization is off, let waLBerla decide which memory allocator
244 * to use. When vectorization is on, the aligned memory allocator is
245 * required, otherwise <tt>cpu_vectorize_info["assume_aligned"]</tt> will
246 * trigger assertions. That is because for single-precision kernels the
247 * waLBerla heuristic in <tt>src/field/allocation/FieldAllocator.h</tt>
248 * will fall back to @c StdFieldAlloc, yet @c AllocateAligned is needed
249 * for intrinsics to work.
250 */
251 template <typename Field>
252 auto add_to_storage(std::string const tag, FloatType value) {
253 auto const &blocks = m_lattice->get_blocks();
254 auto const n_ghost_layers = m_lattice->get_ghost_layers();
255 if constexpr (Architecture == lbmpy::Arch::CPU) {
256 return field::addToStorage<Field>(blocks, tag, FloatType{value},
257 field::fzyx, n_ghost_layers);
258 }
259#if defined(__CUDACC__)
260 else {
261 auto field_id = gpu::addGPUFieldToStorage<GPUField>(
262 blocks, tag, Field::F_SIZE, field::fzyx, n_ghost_layers);
263 if constexpr (std::is_same_v<Field, _DensityField>) {
264 for (auto block = blocks->begin(); block != blocks->end(); ++block) {
265 auto field = block->template getData<GPUField>(field_id);
266 ek::accessor::Scalar::initialize(field, FloatType{value});
267 }
268 } else if constexpr (std::is_same_v<Field, _FluxField>) {
269 for (auto block = blocks->begin(); block != blocks->end(); ++block) {
270 auto field = block->template getData<GPUField>(field_id);
272 std::array<FloatType, FluxCount>{});
273 }
274 }
275 return field_id;
276 }
277#endif
278 }
279
280 void
281 reset_density_boundary_handling(std::shared_ptr<BlockStorage> const &blocks) {
282 auto const [lc, uc] = m_lattice->get_local_grid_range(true);
283 m_boundary_density = std::make_unique<BoundaryModelDensity>(
285 CellInterval{to_cell(lc), to_cell(uc)});
286 }
287
288 void
289 reset_flux_boundary_handling(std::shared_ptr<BlockStorage> const &blocks) {
290 auto const [lc, uc] = m_lattice->get_local_grid_range(true);
291 m_boundary_flux = std::make_shared<BoundaryModelFlux>(
293 CellInterval{to_cell(lc), to_cell(uc)});
294 }
295
297 typename FieldTrait<FloatType, Architecture>::template RegularCommScheme<
298 typename stencil::D3Q27>;
300 typename FieldTrait<FloatType, Architecture>::template BoundaryCommScheme<
301 typename stencil::D3Q27>;
302 std::shared_ptr<FullCommunicator> m_full_communication;
303 std::shared_ptr<BoundaryFullCommunicator> m_boundary_communicator;
304 std::bitset<GhostComm::SIZE> m_pending_ghost_comm;
306 template <class Field>
307 using PackInfo =
309
310public:
311 EKinWalberlaImpl(std::shared_ptr<LatticeWalberla> lattice, double diffusion,
312 double kT, double valency, Utils::Vector3d const &ext_efield,
313 double density, bool advection, bool friction_coupling,
314 bool thermalized, unsigned int seed)
315 : m_diffusion(FloatType_c(diffusion)), m_kT(FloatType_c(kT)),
316 m_valency(FloatType_c(valency)), m_ext_efield(ext_efield),
317 m_advection(advection), m_friction_coupling(friction_coupling),
318 m_seed(seed), m_lattice(std::move(lattice)),
320
321 auto const &blocks = m_lattice->get_blocks();
322 auto const n_ghost_layers = m_lattice->get_ghost_layers();
323
325 add_to_storage<_DensityField>("density field", FloatType_c(density));
327 add_to_storage<_FluxField>("flux field", FloatType_c(0.0));
328
330 std::make_unique<ContinuityKernel>(m_flux_field_id, m_density_field_id);
331
332#if defined(__CUDACC__)
333 m_host_field_allocator =
334 std::make_shared<gpu::HostFieldAllocator<FloatType>>();
335#endif
336
337 if (thermalized) {
338 set_diffusion_kernels(*m_lattice, seed);
339 } else {
340 set_diffusion_kernels();
341 }
342
343 // Init boundary related stuff
344 m_flag_field_density_id = field::addFlagFieldToStorage<FlagField>(
345 blocks, "flag field density", n_ghost_layers);
347
348 m_flag_field_flux_id = field::addFlagFieldToStorage<FlagField>(
349 blocks, "flag field flux", n_ghost_layers);
351
352 m_full_communication = std::make_shared<FullCommunicator>(blocks);
353 m_full_communication->addPackInfo(
354 std::make_shared<PackInfo<DensityField>>(m_density_field_id));
356 std::make_shared<BoundaryFullCommunicator>(blocks);
357 m_boundary_communicator->addPackInfo(
360 auto flux_boundary_packinfo = std::make_shared<
363 flux_boundary_packinfo->setup_boundary_handle(m_lattice, m_boundary_flux);
364 m_boundary_communicator->addPackInfo(flux_boundary_packinfo);
365
367 }
368
369 // Global parameters
370 [[nodiscard]] double get_diffusion() const noexcept override {
371 return m_diffusion;
372 }
373 [[nodiscard]] double get_kT() const noexcept override { return m_kT; }
374 [[nodiscard]] double get_valency() const noexcept override {
375 return m_valency;
376 }
377 [[nodiscard]] bool get_advection() const noexcept override {
378 return m_advection;
379 }
380 [[nodiscard]] bool get_friction_coupling() const noexcept override {
381 return m_friction_coupling;
382 }
383 [[nodiscard]] Utils::Vector3d get_ext_efield() const noexcept override {
384 return m_ext_efield;
385 }
386 [[nodiscard]] bool is_thermalized() const noexcept override {
387 return static_cast<bool>(
388 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux));
389 }
390 [[nodiscard]] unsigned int get_seed() const noexcept override {
391 return m_seed;
392 }
393 [[nodiscard]] std::optional<uint64_t> get_rng_state() const override {
394 auto const kernel =
395 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux);
396 if (!kernel) {
397 return std::nullopt;
398 }
399 return {static_cast<uint64_t>(kernel->getTime_step())};
400 }
401
402 void set_diffusion(double diffusion) override {
403 m_diffusion = FloatType_c(diffusion);
404 auto visitor = [m_diffusion = m_diffusion](auto &kernel) {
405 kernel.setD(m_diffusion);
406 };
407 std::visit(visitor, *m_diffusive_flux);
408 std::visit(visitor, *m_diffusive_flux_electrostatic);
409 }
410
411 void set_kT(double kT) override {
412 m_kT = FloatType_c(kT);
413 std::visit([m_kT = m_kT](auto &kernel) { kernel.setKt(m_kT); },
415 }
416
417 void set_valency(double valency) override {
418 m_valency = FloatType_c(valency);
419 std::visit(
420 [m_valency = m_valency](auto &kernel) { kernel.setZ(m_valency); },
422 }
423
424 void set_advection(bool advection) override { m_advection = advection; }
425
426 void set_friction_coupling(bool friction_coupling) override {
427 m_friction_coupling = friction_coupling;
428 }
429
430 void set_rng_state(uint64_t counter) override {
431 auto const kernel =
432 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux);
433 auto const kernel_electrostatic =
434 std::get_if<DiffusiveFluxKernelElectrostaticThermalized>(
436
437 if (!kernel or !kernel_electrostatic) {
438 throw std::runtime_error("This EK instance is unthermalized");
439 }
440 assert(counter <=
441 static_cast<uint32_t>(std::numeric_limits<uint_t>::max()));
442 kernel->setTime_step(static_cast<uint32_t>(counter));
443 kernel_electrostatic->setTime_step(static_cast<uint32_t>(counter));
444 }
445
446 void set_ext_efield(Utils::Vector3d const &field) override {
447 m_ext_efield = field;
448
449 std::visit(
450 [this](auto &kernel) {
451 kernel.setF_ext_0(FloatType_c(m_ext_efield[0]));
452 kernel.setF_ext_1(FloatType_c(m_ext_efield[1]));
453 kernel.setF_ext_2(FloatType_c(m_ext_efield[2]));
454 },
456 }
457
458 void ghost_communication() override {
461 (*m_full_communication)();
463 }
465 }
466
474
475private:
476 void set_diffusion_kernels() {
477 auto kernel = DiffusiveFluxKernelUnthermalized(
479 m_diffusive_flux = std::make_unique<DiffusiveFluxKernel>(std::move(kernel));
480
481 auto kernel_electrostatic = DiffusiveFluxKernelElectrostaticUnthermalized(
482 m_flux_field_id, BlockDataID{}, m_density_field_id,
483 FloatType_c(m_diffusion), FloatType_c(m_ext_efield[0]),
484 FloatType_c(m_ext_efield[1]), FloatType_c(m_ext_efield[2]),
485 FloatType_c(m_kT), FloatType_c(m_valency));
486
488 std::make_unique<DiffusiveFluxKernelElectrostatic>(
489 std::move(kernel_electrostatic));
490 }
491
492 void set_diffusion_kernels(LatticeWalberla const &lattice,
493 unsigned int seed) {
494 auto const grid_dim = lattice.get_grid_dimensions();
495
496 auto kernel = DiffusiveFluxKernelThermalized(
498 grid_dim[0], grid_dim[1], grid_dim[2], seed, 0);
499
500 auto kernel_electrostatic = DiffusiveFluxKernelElectrostaticThermalized(
501 m_flux_field_id, BlockDataID{}, m_density_field_id,
502 FloatType_c(m_diffusion), FloatType_c(m_ext_efield[0]),
503 FloatType_c(m_ext_efield[1]), FloatType_c(m_ext_efield[2]), grid_dim[0],
504 grid_dim[1], grid_dim[2], FloatType_c(m_kT), seed, 0,
505 FloatType_c(m_valency));
506
507 auto const blocks = lattice.get_blocks();
508
509 for (auto &block : *blocks) {
510 kernel.configure(blocks, &block);
511 kernel_electrostatic.configure(blocks, &block);
512 }
513
514 m_diffusive_flux = std::make_unique<DiffusiveFluxKernel>(std::move(kernel));
516 std::make_unique<DiffusiveFluxKernelElectrostatic>(
517 std::move(kernel_electrostatic));
518 }
519
520 void kernel_boundary_density() {
521 for (auto &block : *m_lattice->get_blocks()) {
522 (*m_boundary_density)(&block);
523 }
524 }
525
526 void kernel_boundary_flux() {
527 for (auto &block : *m_lattice->get_blocks()) {
528 (*m_boundary_flux)(&block);
529 }
530 }
531
532 void kernel_continuity() {
533 for (auto &block : *m_lattice->get_blocks()) {
534 (*m_continuity).run(&block);
535 }
536 }
537
538 void kernel_diffusion() {
539 for (auto &block : *m_lattice->get_blocks()) {
540 std::visit([&block](auto &kernel) { kernel.run(&block); },
542 }
543
544 if (auto *kernel =
545 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux)) {
546 kernel->setTime_step(kernel->getTime_step() + 1u);
547
548 auto *kernel_electrostatic =
549 std::get_if<DiffusiveFluxKernelElectrostaticThermalized>(
551 kernel_electrostatic->setTime_step(kernel_electrostatic->getTime_step() +
552 1u);
553 }
554 }
555
556 void kernel_advection(std::size_t const velocity_id) {
557 auto kernel = AdvectiveFluxKernel(m_flux_field_id, m_density_field_id,
558 BlockDataID(velocity_id));
559 for (auto &block : *m_lattice->get_blocks()) {
560 kernel.run(&block);
561 }
562 }
563
564 void kernel_friction_coupling(std::size_t const force_id,
565 double const lb_density) {
566 auto kernel = FrictionCouplingKernel(
567 BlockDataID(force_id), m_flux_field_id, FloatType_c(get_diffusion()),
568 FloatType_c(get_kT()), FloatType(lb_density));
569 for (auto &block : *m_lattice->get_blocks()) {
570 kernel.run(&block);
571 }
572 }
573
574 void kernel_diffusion_electrostatic(std::size_t const potential_id) {
575 auto const phiID = BlockDataID(potential_id);
576 std::visit([phiID](auto &kernel) { kernel.setPhiID(phiID); },
578
579 for (auto &block : *m_lattice->get_blocks()) {
580 std::visit([&block](auto &kernel) { kernel.run(&block); },
582 }
583
584 if (auto *kernel_electrostatic =
585 std::get_if<DiffusiveFluxKernelElectrostaticThermalized>(
587 kernel_electrostatic->setTime_step(kernel_electrostatic->getTime_step() +
588 1u);
589
590 auto *kernel =
591 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux);
592 kernel->setTime_step(kernel->getTime_step() + 1u);
593 }
594 }
595
596 void kernel_migration() {}
597
598 void update_boundary_fields() {
599 m_boundary_flux->boundary_update();
600 m_boundary_density->boundary_update();
601 }
602
603protected:
604 void integrate_vtk_writers() override {
605 for (auto const &it : m_vtk_auto) {
606 auto &vtk_handle = it.second;
607 if (vtk_handle->enabled) {
608 vtk::writeFiles(vtk_handle->ptr)();
609 vtk_handle->execution_count++;
610 }
611 }
612 }
613
614public:
615 void integrate(std::size_t potential_id, std::size_t velocity_id,
616 std::size_t force_id, double lb_density) override {
617
619 update_boundary_fields();
620
621 if (get_diffusion() == 0.)
622 return;
623
624 if (get_valency() != 0.) {
625 if (potential_id == walberla::BlockDataID{}) {
626 throw std::runtime_error("Walberla EK: electrostatic potential enabled "
627 "but no field accessible. potential id is " +
628 std::to_string(potential_id));
629 }
630 kernel_diffusion_electrostatic(potential_id);
631 } else {
632 kernel_diffusion();
633 }
634
635 kernel_migration();
636 kernel_boundary_flux();
637 // friction coupling
638 if (get_friction_coupling()) {
639 if (force_id == walberla::BlockDataID{}) {
640 throw std::runtime_error("Walberla EK: friction coupling enabled but "
641 "no force field accessible. force_id is " +
642 std::to_string(force_id) +
643 ". Hint: LB may be inactive.");
644 }
645 kernel_friction_coupling(force_id, lb_density);
646 }
647
648 if (get_advection()) {
649 if (velocity_id == walberla::BlockDataID{}) {
650 throw std::runtime_error("Walberla EK: advection enabled but no "
651 "velocity field accessible. velocity_id is " +
652 std::to_string(velocity_id) +
653 ". Hint: LB may be inactive.");
654 }
655 kernel_advection(velocity_id);
656 kernel_boundary_flux();
657 }
658 kernel_continuity();
659
660 // is this the expected behavior when reactions are included?
661 kernel_boundary_density();
663
664 // Handle VTK writers
666 }
667
668 [[nodiscard]] std::size_t get_density_id() const noexcept override {
669 static_assert(std::is_same_v<std::size_t, walberla::uint_t>);
670 return static_cast<std::size_t>(m_density_field_id);
671 }
672
673 bool set_node_density(Utils::Vector3i const &node, double density) override {
675 auto bc = get_block_and_cell(get_lattice(), node, false);
676 if (!bc)
677 return false;
678
679 auto density_field =
680 bc->block->template getData<DensityField>(m_density_field_id);
681 ek::accessor::Scalar::set(density_field, FloatType_c(density), bc->cell);
682 return true;
683 }
684
685 [[nodiscard]] std::optional<double>
687 bool consider_ghosts = false) const override {
688 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
689
690 if (!bc)
691 return std::nullopt;
692
693 auto const density_field =
694 bc->block->template getData<DensityField>(m_density_field_id);
695 return {double_c(ek::accessor::Scalar::get(density_field, bc->cell))};
696 }
697
698 [[nodiscard]] std::vector<double>
700 Utils::Vector3i const &upper_corner) const override {
701 std::vector<double> out;
702 uint_t values_size = 0;
703 auto const &lattice = get_lattice();
704 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
705 out = std::vector<double>(ci->numCells());
706 for (auto &block : *lattice.get_blocks()) {
707 auto const block_offset = lattice.get_block_corner(block, true);
708 if (auto const bci = get_block_interval(
709 lattice, lower_corner, upper_corner, block_offset, block)) {
710 auto const density_field =
711 block.template getData<DensityField>(m_density_field_id);
712 auto const values = ek::accessor::Scalar::get(density_field, *bci);
713 assert(values.size() == bci->numCells());
714 values_size += bci->numCells();
715 auto kernel = [&values, &out](unsigned const block_index,
716 unsigned const local_index,
717 Utils::Vector3i const &) {
718 out[local_index] = double_c(values[block_index]);
719 };
720
721 copy_block_buffer(*bci, *ci, block_offset, lower_corner, kernel);
722 }
723 }
724 assert(values_size == ci->numCells());
725 }
726 return out;
727 }
728
729 void set_slice_density(Utils::Vector3i const &lower_corner,
730 Utils::Vector3i const &upper_corner,
731 std::vector<double> const &density) override {
733 auto const &lattice = get_lattice();
734 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
735 assert(density.size() == ci->numCells());
736 for (auto &block : *lattice.get_blocks()) {
737 auto const block_offset = lattice.get_block_corner(block, true);
738 if (auto const bci = get_block_interval(
739 lattice, lower_corner, upper_corner, block_offset, block)) {
740 auto const density_field =
741 block.template getData<DensityField>(m_density_field_id);
742 std::vector<FloatType> values(bci->numCells());
743
744 auto kernel = [&values, &density](unsigned const block_index,
745 unsigned const local_index,
746 Utils::Vector3i const &) {
747 values[block_index] = numeric_cast<FloatType>(density[local_index]);
748 };
749
750 copy_block_buffer(*bci, *ci, block_offset, lower_corner, kernel);
751 ek::accessor::Scalar::set(density_field, values, *bci);
752 }
753 }
754 }
755 }
756
757 [[nodiscard]] std::optional<Utils::Vector3d>
759 bool consider_ghosts = false) const override {
760 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
761
762 if (!bc)
763 return std::nullopt;
764
765 auto const flux_field =
766 bc->block->template getData<FluxField>(m_flux_field_id);
767 return to_vector3d(ek::accessor::Flux::get_vector(flux_field, bc->cell));
768 }
769
770 std::vector<double>
772 Utils::Vector3i const &upper_corner) const override {
773 std::vector<double> out;
774 uint_t values_size = 0;
775 auto const &lattice = get_lattice();
776 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
777 out = std::vector<double>(3u * ci->numCells());
778 for (auto &block : *lattice.get_blocks()) {
779 auto const block_offset = lattice.get_block_corner(block, true);
780 if (auto const bci = get_block_interval(
781 lattice, lower_corner, upper_corner, block_offset, block)) {
782 auto const flux_field =
783 block.template getData<FluxField>(m_flux_field_id);
784 auto const values = ek::accessor::Flux::get_vector(flux_field, *bci);
785 assert(values.size() == 3u * bci->numCells());
786 values_size += 3u * bci->numCells();
787
788 auto kernel = [&values, &out, this](unsigned const block_index,
789 unsigned const local_index,
790 Utils::Vector3i const &node) {
791 if (m_boundary_flux->node_is_boundary(node)) {
792 auto const &vec =
793 m_boundary_flux->get_node_value_at_boundary(node);
794 for (uint_t f = 0u; f < 3u; ++f) {
795 out[3u * local_index + f] = double_c(vec[f]);
796 }
797 } else {
798 for (uint_t f = 0u; f < 3u; ++f) {
799 out[3u * local_index + f] =
800 double_c(values[3u * block_index + f]);
801 }
802 }
803 };
804
805 copy_block_buffer(*bci, *ci, block_offset, lower_corner, kernel);
806 }
807 }
808 assert(values_size == 3u * ci->numCells());
809 }
810 return out;
811 }
812
817
818 void clear_density_boundaries() override {
820 }
821
823 Utils::Vector3d const &flux) override {
825 auto bc = get_block_and_cell(get_lattice(), node, true);
826 if (!bc)
827 return false;
828
829 m_boundary_flux->set_node_value_at_boundary(
830 node, to_vector3<FloatType>(flux), *bc);
831 return true;
832 }
833
834 [[nodiscard]] std::optional<Utils::Vector3d>
836 bool consider_ghosts = false) const override {
837 assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::FLB)));
838 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
839 if (!bc or !m_boundary_flux->node_is_boundary(node))
840 return std::nullopt;
841
842 return {to_vector3d(m_boundary_flux->get_node_value_at_boundary(node))};
843 }
844
847 auto bc = get_block_and_cell(get_lattice(), node, true);
848 if (!bc)
849 return false;
850
851 m_boundary_flux->remove_node_from_boundary(node, *bc);
852 return true;
853 }
854
856 double density) override {
857 auto bc = get_block_and_cell(get_lattice(), node, true);
858 if (!bc)
859 return false;
860
861 m_boundary_density->set_node_value_at_boundary(node, FloatType_c(density),
862 *bc);
863
864 return true;
865 }
866
867 [[nodiscard]] std::optional<double>
869 bool consider_ghosts = false) const override {
870 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
871 if (!bc or !m_boundary_density->node_is_boundary(node))
872 return std::nullopt;
873
874 return {double_c(m_boundary_density->get_node_value_at_boundary(node))};
875 }
876
878 Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
879 std::vector<std::optional<double>> const &density) override {
880 auto const &lattice = get_lattice();
881 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
882 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
883 auto const lower_cell = ci->min();
884 auto const upper_cell = ci->max();
885 auto it = density.begin();
886 assert(density.size() == ci->numCells());
887 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
888 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
889 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
890 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
891 auto const bc = get_block_and_cell(lattice, node, false);
892 auto const &opt = *it;
893 if (opt) {
894 m_boundary_density->set_node_value_at_boundary(
895 node, FloatType_c(*opt), *bc);
896 } else {
897 m_boundary_density->remove_node_from_boundary(node, *bc);
898 }
899 ++it;
900 }
901 }
902 }
903 }
904 }
905
906 [[nodiscard]] std::vector<std::optional<double>>
908 Utils::Vector3i const &lower_corner,
909 Utils::Vector3i const &upper_corner) const override {
910 std::vector<std::optional<double>> out;
911 auto const &lattice = get_lattice();
912 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
913 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
914 auto const lower_cell = ci->min();
915 auto const upper_cell = ci->max();
916 auto const n_values = ci->numCells();
917 out.reserve(n_values);
918 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
919 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
920 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
921 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
922 if (m_boundary_density->node_is_boundary(node)) {
923 out.emplace_back(double_c(
924 m_boundary_density->get_node_value_at_boundary(node)));
925 } else {
926 out.emplace_back(std::nullopt);
927 }
928 }
929 }
930 }
931 assert(out.size() == n_values);
932 }
933 return out;
934 }
935
937 Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
938 std::vector<std::optional<Utils::Vector3d>> const &flux) override {
940 auto const &lattice = get_lattice();
941 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
942 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
943 auto const lower_cell = ci->min();
944 auto const upper_cell = ci->max();
945 auto it = flux.begin();
946 assert(flux.size() == ci->numCells());
947 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
948 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
949 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
950 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
951 auto const bc = get_block_and_cell(lattice, node, false);
952 auto const &opt = *it;
953 if (opt) {
954 m_boundary_flux->set_node_value_at_boundary(
955 node, to_vector3<FloatType>(*opt), *bc);
956 } else {
957 m_boundary_flux->remove_node_from_boundary(node, *bc);
958 }
959 ++it;
960 }
961 }
962 }
963 }
964 }
965
966 [[nodiscard]] std::vector<std::optional<Utils::Vector3d>>
968 Utils::Vector3i const &lower_corner,
969 Utils::Vector3i const &upper_corner) const override {
970 std::vector<std::optional<Utils::Vector3d>> out;
971 auto const &lattice = get_lattice();
972 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
973 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
974 auto const lower_cell = ci->min();
975 auto const upper_cell = ci->max();
976 auto const n_values = ci->numCells();
977 out.reserve(n_values);
978 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
979 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
980 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
981 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
982 if (m_boundary_flux->node_is_boundary(node)) {
983 out.emplace_back(to_vector3d(
984 m_boundary_flux->get_node_value_at_boundary(node)));
985 } else {
986 out.emplace_back(std::nullopt);
987 }
988 }
989 }
990 }
991 assert(out.size() == n_values);
992 }
993 return out;
994 }
995
996 [[nodiscard]] std::vector<bool>
998 Utils::Vector3i const &upper_corner) const override {
999 std::vector<bool> out;
1000 auto const &lattice = get_lattice();
1001 if (auto const ci = get_interval(lattice, lower_corner, upper_corner)) {
1002 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
1003 auto const lower_cell = ci->min();
1004 auto const upper_cell = ci->max();
1005 auto const n_values = ci->numCells();
1006 out.reserve(n_values);
1007 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
1008 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
1009 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
1010 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
1011 out.emplace_back(m_boundary_density->node_is_boundary(node) or
1012 m_boundary_flux->node_is_boundary(node));
1013 }
1014 }
1015 }
1016 assert(out.size() == n_values);
1017 }
1018 return out;
1019 }
1020
1022 auto bc = get_block_and_cell(get_lattice(), node, true);
1023 if (!bc)
1024 return false;
1025
1026 m_boundary_density->remove_node_from_boundary(node, *bc);
1027
1028 return true;
1029 }
1030
1031 [[nodiscard]] std::optional<bool>
1033 bool consider_ghosts) const override {
1034 assert(not(consider_ghosts and m_pending_ghost_comm.test(GhostComm::FLB)));
1035 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
1036 if (!bc)
1037 return std::nullopt;
1038
1039 return {m_boundary_flux->node_is_boundary(node)};
1040 }
1041
1042 [[nodiscard]] std::optional<bool>
1044 bool consider_ghosts) const override {
1045 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
1046 if (!bc)
1047 return std::nullopt;
1048
1049 return {m_boundary_density->node_is_boundary(node)};
1050 }
1051
1052 [[nodiscard]] std::optional<bool>
1054 bool consider_ghosts = false) const override {
1055 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
1056 if (!bc)
1057 return std::nullopt;
1058
1059 return {m_boundary_density->node_is_boundary(node) or
1060 m_boundary_flux->node_is_boundary(node)};
1061 }
1062
1064 const std::vector<int> &raster_flat,
1065 const std::vector<double> &data_flat) override {
1067 auto const grid_size = get_lattice().get_grid_dimensions();
1068 auto const data = fill_3D_vector_array(data_flat, grid_size);
1069 set_boundary_from_grid(*m_boundary_flux, get_lattice(), raster_flat, data);
1071 }
1072
1074 const std::vector<int> &raster_flat,
1075 const std::vector<double> &data_flat) override {
1076 auto const grid_size = get_lattice().get_grid_dimensions();
1077 auto const data = fill_3D_scalar_array(data_flat, grid_size);
1079 data);
1081 }
1082
1084
1086 m_boundary_density->boundary_update();
1087 }
1088
1089 [[nodiscard]] LatticeWalberla const &get_lattice() const noexcept override {
1090 return *m_lattice;
1091 }
1092
1093 [[nodiscard]] bool is_gpu() const noexcept override {
1094 return Architecture == lbmpy::Arch::GPU;
1095 }
1096
1097 void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override {
1098 field::FlagFieldCellFilter<FlagField> fluid_filter(m_flag_field_density_id);
1099 fluid_filter.addFlag(Boundary_flag);
1100 vtk_obj.addCellExclusionFilter(fluid_filter);
1101 }
1102
1103protected:
1104 template <typename Field_T, uint_t F_SIZE_ARG, typename OutputType>
1105 class VTKWriter : public vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG> {
1106 public:
1107 VTKWriter(ConstBlockDataID const &block_id, std::string const &id,
1108 FloatType unit_conversion)
1109 : vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG>(id),
1110 m_block_id(block_id), m_field(nullptr),
1111 m_conversion(unit_conversion) {}
1112
1113 protected:
1114 void configure() override {
1115 WALBERLA_ASSERT_NOT_NULLPTR(this->block_);
1116 m_field = this->block_->template getData<Field_T>(m_block_id);
1117 }
1118
1119 ConstBlockDataID const m_block_id;
1120 Field_T const *m_field;
1121 FloatType const m_conversion;
1122 };
1123
1124#if defined(__CUDACC__)
1125 template <typename OutputType = float>
1126 class DensityVTKWriter : public VTKWriter<DensityFieldCpu, 1u, OutputType> {
1127 public:
1129 using Base::Base;
1130 using Base::evaluate;
1131
1132 protected:
1133 OutputType evaluate(cell_idx_t const x, cell_idx_t const y,
1134 cell_idx_t const z, cell_idx_t const) override {
1135 WALBERLA_ASSERT_NOT_NULLPTR(this->m_field);
1136 auto const density = ek::accessor::Scalar::get(this->m_field, {x, y, z});
1137 return numeric_cast<OutputType>(this->m_conversion * density);
1138 }
1139 };
1140#else
1141 template <typename OutputType = float>
1142 class DensityVTKWriter : public VTKWriter<DensityField, 1u, OutputType> {
1143 public:
1145 using Base::Base;
1146 using Base::evaluate;
1147
1148 protected:
1149 OutputType evaluate(cell_idx_t const x, cell_idx_t const y,
1150 cell_idx_t const z, cell_idx_t const) override {
1151 WALBERLA_ASSERT_NOT_NULLPTR(this->m_field);
1152 auto const density = ek::accessor::Scalar::get(this->m_field, {x, y, z});
1153 return numeric_cast<OutputType>(this->m_conversion * density);
1154 }
1155 };
1156#endif
1157
1158#if defined(__CUDACC__)
1159 template <typename OutputType = float>
1160 class FluxVTKWriter : public VTKWriter<FluxFieldCpu, 3u, OutputType> {
1161 public:
1162 using Base = VTKWriter<FluxFieldCpu, 3u, OutputType>;
1163 using Base::Base;
1164 using Base::evaluate;
1165
1166 protected:
1167 OutputType evaluate(cell_idx_t const x, cell_idx_t const y,
1168 cell_idx_t const z, cell_idx_t const f) override {
1169 WALBERLA_ASSERT_NOT_NULLPTR(this->m_field);
1170 auto const flux =
1171 ek::accessor::Flux::get_vector(this->m_field, {x, y, z});
1172 return numeric_cast<OutputType>(this->m_conversion * flux[uint_c(f)]);
1173 }
1174 };
1175#else
1176 template <typename OutputType = float>
1177 class FluxVTKWriter : public VTKWriter<FluxField, 3u, OutputType> {
1178 public:
1180 using Base::Base;
1181 using Base::evaluate;
1182
1183 protected:
1184 OutputType evaluate(cell_idx_t const x, cell_idx_t const y,
1185 cell_idx_t const z, cell_idx_t const f) override {
1186 WALBERLA_ASSERT_NOT_NULLPTR(this->m_field);
1187 auto const flux =
1188 ek::accessor::Flux::get_vector(this->m_field, {x, y, z});
1189 return numeric_cast<OutputType>(this->m_conversion * flux[uint_c(f)]);
1190 }
1191 };
1192#endif
1193
1194public:
1195 void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj,
1196 LatticeModel::units_map const &units,
1197 int flag_observables) override {
1198#if defined(__CUDACC__)
1199 auto const allocate_cpu_field_if_empty =
1200 [&]<typename Field>(auto const &blocks, std::string name,
1201 std::optional<BlockDataID> &cpu_field) {
1202 if (not cpu_field) {
1203 cpu_field = field::addToStorage<Field>(
1204 blocks, name, FloatType{0}, field::fzyx,
1205 m_lattice->get_ghost_layers(), m_host_field_allocator);
1206 }
1207 };
1208#endif
1209 if (flag_observables & static_cast<int>(EKOutputVTK::density)) {
1210 auto const unit_conversion = FloatType_c(units.at("density"));
1211#if defined(__CUDACC__)
1212 if constexpr (Architecture == lbmpy::Arch::GPU) {
1213 auto const &blocks = m_lattice->get_blocks();
1214 allocate_cpu_field_if_empty.template operator()<DensityFieldCpu>(
1215 blocks, "density_cpu", m_density_cpu_field_id);
1216 vtk_obj.addBeforeFunction(
1217 gpu::fieldCpyFunctor<DensityFieldCpu, DensityField>(
1218 blocks, *m_density_cpu_field_id, m_density_field_id));
1219 vtk_obj.addCellDataWriter(make_shared<DensityVTKWriter<float>>(
1220 *m_density_cpu_field_id, "density", unit_conversion));
1221 } else {
1222#endif
1223 vtk_obj.addCellDataWriter(make_shared<DensityVTKWriter<float>>(
1224 m_density_field_id, "density", unit_conversion));
1225#if defined(__CUDACC__)
1226 }
1227#endif
1228 }
1229 if (flag_observables & static_cast<int>(EKOutputVTK::flux)) {
1230 auto const unit_conversion = FloatType_c(units.at("flux"));
1231#if defined(__CUDACC__)
1232 if constexpr (Architecture == lbmpy::Arch::GPU) {
1233 auto const &blocks = m_lattice->get_blocks();
1234 allocate_cpu_field_if_empty.template operator()<FluxFieldCpu>(
1235 blocks, "flux_cpu", m_flux_cpu_field_id);
1236 vtk_obj.addBeforeFunction(gpu::fieldCpyFunctor<FluxFieldCpu, FluxField>(
1237 blocks, *m_flux_cpu_field_id, m_flux_field_id));
1238 vtk_obj.addCellDataWriter(make_shared<FluxVTKWriter<float>>(
1239 *m_flux_cpu_field_id, "flux", unit_conversion));
1240 } else {
1241#endif
1242 vtk_obj.addCellDataWriter(make_shared<FluxVTKWriter<float>>(
1243 m_flux_field_id, "flux", unit_conversion));
1244#if defined(__CUDACC__)
1245 }
1246#endif
1247 }
1248 }
1249
1250 ~EKinWalberlaImpl() override = default;
1251};
1252
1253} // 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.
VTKWriter< DensityField, 1u, OutputType > Base
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
VTKWriter< FluxField, 3u, OutputType > Base
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:176
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
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)
void copy_block_buffer(CellInterval const &bci, CellInterval const &ci, Utils::Vector3i const &block_offset, Utils::Vector3i const &lower_corner, Kernel &&kernel)
Synchronize data between a sliced block and a container.
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.