Loading [MathJax]/extensions/TeX/AMSmath.js
ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages Concepts
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
32#include "../BoundaryHandling.hpp"
33#include "../utils/boundary.hpp"
34#include "../utils/types_conversion.hpp"
35#include "ek_kernels.hpp"
36
40
41#include <utils/Vector.hpp>
42
43#include <cstddef>
44#include <iterator>
45#include <memory>
46#include <optional>
47#include <stdexcept>
48#include <string>
49#include <type_traits>
50#include <variant>
51#include <vector>
52
53namespace walberla {
54
55/** @brief Class that runs and controls the EK on waLBerla. */
56template <std::size_t FluxCount = 13, typename FloatType = double>
58 using ContinuityKernel =
60 using DiffusiveFluxKernelUnthermalized =
62 using DiffusiveFluxKernelThermalized =
64 using AdvectiveFluxKernel =
66 using FrictionCouplingKernel =
68 using DiffusiveFluxKernelElectrostaticUnthermalized =
70 using DiffusiveFluxKernelElectrostaticThermalized =
71 typename detail::KernelTrait<
72 FloatType>::DiffusiveFluxKernelElectrostaticThermalized;
73
74 using DiffusiveFluxKernel = std::variant<DiffusiveFluxKernelUnthermalized,
75 DiffusiveFluxKernelThermalized>;
76 using DiffusiveFluxKernelElectrostatic =
77 std::variant<DiffusiveFluxKernelElectrostaticUnthermalized,
78 DiffusiveFluxKernelElectrostaticThermalized>;
79
80 using Dirichlet = typename detail::KernelTrait<FloatType>::Dirichlet;
81 using FixedFlux = typename detail::KernelTrait<FloatType>::FixedFlux;
82
83protected:
84 // Type definitions
85 using FluxField = GhostLayerField<FloatType, FluxCount>;
86 using FlagField = walberla::FlagField<walberla::uint8_t>;
87 using DensityField = GhostLayerField<FloatType, 1>;
88
91
93
94public:
95 template <typename T> FloatType FloatType_c(T t) {
96 return numeric_cast<FloatType>(t);
97 }
98
99 [[nodiscard]] std::size_t stencil_size() const noexcept override {
100 return FluxCount;
101 }
102
103 [[nodiscard]] bool is_double_precision() const noexcept override {
104 return std::is_same<FloatType, double>::value;
105 }
106
107private:
108 FloatType m_diffusion;
109 FloatType m_kT;
110 FloatType m_valency;
111 Utils::Vector3d m_ext_efield;
112 bool m_advection;
113 bool m_friction_coupling;
114 unsigned int m_seed;
115
116protected:
117 // Block data access handles
120
121 BlockDataID m_flux_field_id;
123
126
127 /** Flag for domain cells, i.e. all cells. */
128 FlagUID const Domain_flag{"domain"};
129 /** Flag for boundary cells. */
130 FlagUID const Boundary_flag{"boundary"};
131
132 /** Block forest */
133 std::shared_ptr<LatticeWalberla> m_lattice;
134
135 std::unique_ptr<BoundaryModelDensity> m_boundary_density;
136 std::unique_ptr<BoundaryModelFlux> m_boundary_flux;
137
138 std::unique_ptr<DiffusiveFluxKernel> m_diffusive_flux;
139 std::unique_ptr<DiffusiveFluxKernelElectrostatic>
141 std::unique_ptr<ContinuityKernel> m_continuity;
142
143 // ResetFlux + external force
144 // TODO: kernel for that
145 // std::shared_ptr<ResetForce<PdfField, VectorField>> m_reset_force;
146
147 [[nodiscard]] std::optional<CellInterval>
148 get_interval(Utils::Vector3i const &lower_corner,
149 Utils::Vector3i const &upper_corner) const {
150 auto const &lattice = get_lattice();
151 auto const &cell_min = lower_corner;
152 auto const cell_max = upper_corner - Utils::Vector3i::broadcast(1);
153 auto const lower_bc = get_block_and_cell(lattice, cell_min, true);
154 auto const upper_bc = get_block_and_cell(lattice, cell_max, true);
155 if (not lower_bc or not upper_bc) {
156 return std::nullopt;
157 }
158 assert(&(*(lower_bc->block)) == &(*(upper_bc->block)));
159 return {CellInterval(lower_bc->cell, upper_bc->cell)};
160 }
161
162 void
163 reset_density_boundary_handling(std::shared_ptr<BlockStorage> const &blocks) {
164 m_boundary_density = std::make_unique<BoundaryModelDensity>(
166 }
167
168 void
169 reset_flux_boundary_handling(std::shared_ptr<BlockStorage> const &blocks) {
170 m_boundary_flux = std::make_unique<BoundaryModelFlux>(
172 }
173
174 using FullCommunicator = blockforest::communication::UniformBufferedScheme<
175 typename stencil::D3Q27>;
176 std::shared_ptr<FullCommunicator> m_full_communication;
177
178public:
179 EKinWalberlaImpl(std::shared_ptr<LatticeWalberla> lattice, double diffusion,
180 double kT, double valency, Utils::Vector3d const &ext_efield,
181 double density, bool advection, bool friction_coupling,
182 bool thermalized, unsigned int seed)
183 : m_diffusion(FloatType_c(diffusion)), m_kT(FloatType_c(kT)),
184 m_valency(FloatType_c(valency)), m_ext_efield(ext_efield),
185 m_advection(advection), m_friction_coupling(friction_coupling),
186 m_seed(seed), m_lattice(std::move(lattice)) {
187
188 auto const &blocks = m_lattice->get_blocks();
189 auto const n_ghost_layers = m_lattice->get_ghost_layers();
190
191 m_density_field_id = field::addToStorage<DensityField>(
192 blocks, "density field", FloatType_c(density), field::fzyx,
193 n_ghost_layers);
195 field::addFlattenedShallowCopyToStorage<DensityField>(
196 blocks, m_density_field_id, "flattened density field");
197 m_flux_field_id = field::addToStorage<FluxField>(
198 blocks, "flux field", FloatType{0}, field::fzyx, n_ghost_layers);
200 field::addFlattenedShallowCopyToStorage<FluxField>(
201 blocks, m_flux_field_id, "flattened flux field");
202
203 m_continuity = std::make_unique<ContinuityKernel>(
205
206 if (thermalized) {
207 set_diffusion_kernels(*m_lattice, seed);
208 } else {
209 set_diffusion_kernels();
210 }
211
212 // Init boundary related stuff
213 m_flag_field_density_id = field::addFlagFieldToStorage<FlagField>(
214 blocks, "flag field density", n_ghost_layers);
216
217 m_flag_field_flux_id = field::addFlagFieldToStorage<FlagField>(
218 blocks, "flag field flux", n_ghost_layers);
220
221 m_full_communication = std::make_shared<FullCommunicator>(blocks);
222 m_full_communication->addPackInfo(
223 std::make_shared<field::communication::PackInfo<DensityField>>(
225 }
226
227 // Global parameters
228 [[nodiscard]] double get_diffusion() const noexcept override {
229 return m_diffusion;
230 }
231 [[nodiscard]] double get_kT() const noexcept override { return m_kT; }
232 [[nodiscard]] double get_valency() const noexcept override {
233 return m_valency;
234 }
235 [[nodiscard]] bool get_advection() const noexcept override {
236 return m_advection;
237 }
238 [[nodiscard]] bool get_friction_coupling() const noexcept override {
239 return m_friction_coupling;
240 }
241 [[nodiscard]] Utils::Vector3d get_ext_efield() const noexcept override {
242 return m_ext_efield;
243 }
244 [[nodiscard]] bool is_thermalized() const noexcept override {
245 return static_cast<bool>(
246 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux));
247 }
248 [[nodiscard]] unsigned int get_seed() const noexcept override {
249 return m_seed;
250 }
251 [[nodiscard]] std::optional<uint64_t> get_rng_state() const override {
252 auto const kernel =
253 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux);
254 if (!kernel) {
255 return std::nullopt;
256 }
257 return {static_cast<uint64_t>(kernel->getTime_step())};
258 }
259
260 void set_diffusion(double diffusion) override {
261 m_diffusion = FloatType_c(diffusion);
262 auto visitor = [m_diffusion = m_diffusion](auto &kernel) {
263 kernel.setD(m_diffusion);
264 };
265 std::visit(visitor, *m_diffusive_flux);
266 std::visit(visitor, *m_diffusive_flux_electrostatic);
267 }
268
269 void set_kT(double kT) override {
270 m_kT = FloatType_c(kT);
271 std::visit([m_kT = m_kT](auto &kernel) { kernel.setKt(m_kT); },
273 }
274
275 void set_valency(double valency) override {
276 m_valency = FloatType_c(valency);
277 std::visit(
278 [m_valency = m_valency](auto &kernel) { kernel.setZ(m_valency); },
280 }
281
282 void set_advection(bool advection) override { m_advection = advection; }
283
284 void set_friction_coupling(bool friction_coupling) override {
285 m_friction_coupling = friction_coupling;
286 }
287
288 void set_rng_state(uint64_t counter) override {
289 auto const kernel =
290 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux);
291 auto const kernel_electrostatic =
292 std::get_if<DiffusiveFluxKernelElectrostaticThermalized>(
294
295 if (!kernel or !kernel_electrostatic) {
296 throw std::runtime_error("This EK instance is unthermalized");
297 }
298 assert(counter <=
299 static_cast<uint32_t>(std::numeric_limits<uint_t>::max()));
300 kernel->setTime_step(static_cast<uint32_t>(counter));
301 kernel_electrostatic->setTime_step(static_cast<uint32_t>(counter));
302 }
303
304 void set_ext_efield(Utils::Vector3d const &field) override {
305 m_ext_efield = field;
306
307 std::visit(
308 [this](auto &kernel) {
309 kernel.setF_ext_0(FloatType_c(m_ext_efield[0]));
310 kernel.setF_ext_1(FloatType_c(m_ext_efield[1]));
311 kernel.setF_ext_2(FloatType_c(m_ext_efield[2]));
312 },
314 }
315
316 void ghost_communication() override { (*m_full_communication)(); }
317
318private:
319 void set_diffusion_kernels() {
320 auto kernel = DiffusiveFluxKernelUnthermalized(m_flux_field_flattened_id,
322 FloatType_c(m_diffusion));
323 m_diffusive_flux = std::make_unique<DiffusiveFluxKernel>(std::move(kernel));
324
325 auto kernel_electrostatic = DiffusiveFluxKernelElectrostaticUnthermalized(
327 FloatType_c(m_diffusion), FloatType_c(m_ext_efield[0]),
328 FloatType_c(m_ext_efield[1]), FloatType_c(m_ext_efield[2]),
329 FloatType_c(m_kT), FloatType_c(m_valency));
330
332 std::make_unique<DiffusiveFluxKernelElectrostatic>(
333 std::move(kernel_electrostatic));
334 }
335
336 void set_diffusion_kernels(LatticeWalberla const &lattice,
337 unsigned int seed) {
338 auto const grid_dim = lattice.get_grid_dimensions();
339
340 auto kernel = DiffusiveFluxKernelThermalized(
342 FloatType_c(m_diffusion), grid_dim[0], grid_dim[1], grid_dim[2], seed,
343 0);
344
345 auto kernel_electrostatic = DiffusiveFluxKernelElectrostaticThermalized(
347 FloatType_c(m_diffusion), FloatType_c(m_ext_efield[0]),
348 FloatType_c(m_ext_efield[1]), FloatType_c(m_ext_efield[2]), grid_dim[0],
349 grid_dim[1], grid_dim[2], FloatType_c(m_kT), seed, 0,
350 FloatType_c(m_valency));
351
352 auto const blocks = lattice.get_blocks();
353
354 for (auto &block : *blocks) {
355 kernel.configure(blocks, &block);
356 kernel_electrostatic.configure(blocks, &block);
357 }
358
359 m_diffusive_flux = std::make_unique<DiffusiveFluxKernel>(std::move(kernel));
361 std::make_unique<DiffusiveFluxKernelElectrostatic>(
362 std::move(kernel_electrostatic));
363 }
364
365 void kernel_boundary_density() {
366 for (auto &block : *m_lattice->get_blocks()) {
367 (*m_boundary_density)(&block);
368 }
369 }
370
371 void kernel_boundary_flux() {
372 for (auto &block : *m_lattice->get_blocks()) {
373 (*m_boundary_flux)(&block);
374 }
375 }
376
377 void kernel_continuity() {
378 for (auto &block : *m_lattice->get_blocks()) {
379 (*m_continuity).run(&block);
380 }
381 }
382
383 void kernel_diffusion() {
384 for (auto &block : *m_lattice->get_blocks()) {
385 std::visit([&block](auto &kernel) { kernel.run(&block); },
387 }
388
389 if (auto *kernel =
390 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux)) {
391 kernel->setTime_step(kernel->getTime_step() + 1u);
392
393 auto *kernel_electrostatic =
394 std::get_if<DiffusiveFluxKernelElectrostaticThermalized>(
396 kernel_electrostatic->setTime_step(kernel_electrostatic->getTime_step() +
397 1u);
398 }
399 }
400
401 void kernel_advection(const std::size_t &velocity_id) {
402 auto kernel =
404 BlockDataID(velocity_id));
405 for (auto &block : *m_lattice->get_blocks()) {
406 kernel.run(&block);
407 }
408 }
409
410 void kernel_friction_coupling(const std::size_t &force_id) {
411 auto kernel = FrictionCouplingKernel(
412 BlockDataID(force_id), m_flux_field_flattened_id,
414 for (auto &block : *m_lattice->get_blocks()) {
415 kernel.run(&block);
416 }
417 }
418
419 void kernel_diffusion_electrostatic(const std::size_t &potential_id) {
420 auto const phiID = BlockDataID(potential_id);
421 std::visit([phiID](auto &kernel) { kernel.setPhiID(phiID); },
423
424 for (auto &block : *m_lattice->get_blocks()) {
425 std::visit([&block](auto &kernel) { kernel.run(&block); },
427 }
428
429 if (auto *kernel_electrostatic =
430 std::get_if<DiffusiveFluxKernelElectrostaticThermalized>(
432 kernel_electrostatic->setTime_step(kernel_electrostatic->getTime_step() +
433 1u);
434
435 auto *kernel =
436 std::get_if<DiffusiveFluxKernelThermalized>(&*m_diffusive_flux);
437 kernel->setTime_step(kernel->getTime_step() + 1u);
438 }
439 }
440
441 void kernel_migration() {}
442
443 void updated_boundary_fields() {
444 m_boundary_flux->boundary_update();
445 m_boundary_density->boundary_update();
446 }
447
448protected:
449 void integrate_vtk_writers() override {
450 for (auto const &it : m_vtk_auto) {
451 auto &vtk_handle = it.second;
452 if (vtk_handle->enabled) {
453 vtk::writeFiles(vtk_handle->ptr)();
454 vtk_handle->execution_count++;
455 }
456 }
457 }
458
459public:
460 void integrate(std::size_t potential_id, std::size_t velocity_id,
461 std::size_t force_id) override {
462
463 updated_boundary_fields();
464
465 if (get_diffusion() == 0.)
466 return;
467
468 if (get_valency() != 0.) {
469 if (potential_id == walberla::BlockDataID{}) {
470 throw std::runtime_error("Walberla EK: electrostatic potential enabled "
471 "but no field accessible. potential id is " +
472 std::to_string(potential_id));
473 }
474 kernel_diffusion_electrostatic(potential_id);
475 } else {
476 kernel_diffusion();
477 }
478
479 kernel_migration();
480 kernel_boundary_flux();
481 // friction coupling
482 if (get_friction_coupling()) {
483 if (force_id == walberla::BlockDataID{}) {
484 throw std::runtime_error("Walberla EK: friction coupling enabled but "
485 "no force field accessible. force_id is " +
486 std::to_string(force_id) +
487 ". Hint: LB may be inactive.");
488 }
489 kernel_friction_coupling(force_id);
490 }
491
492 if (get_advection()) {
493 if (velocity_id == walberla::BlockDataID{}) {
494 throw std::runtime_error("Walberla EK: advection enabled but no "
495 "velocity field accessible. velocity_id is " +
496 std::to_string(velocity_id) +
497 ". Hint: LB may be inactive.");
498 }
499 kernel_advection(velocity_id);
500 }
501 kernel_continuity();
502
503 // is this the expected behavior when reactions are included?
504 kernel_boundary_density();
505
506 // Handle VTK writers
508 }
509
510 [[nodiscard]] std::size_t get_density_id() const noexcept override {
511 static_assert(std::is_same_v<std::size_t, walberla::uint_t>);
512 return static_cast<std::size_t>(m_density_field_id);
513 }
514
515 bool set_node_density(Utils::Vector3i const &node, double density) override {
516 auto bc = get_block_and_cell(get_lattice(), node, false);
517 if (!bc)
518 return false;
519
520 auto density_field =
521 bc->block->template getData<DensityField>(m_density_field_id);
522 density_field->get(bc->cell) = FloatType_c(density);
523
524 return true;
525 }
526
527 [[nodiscard]] std::optional<double>
529 bool consider_ghosts = false) const override {
530 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
531
532 if (!bc)
533 return std::nullopt;
534
535 auto const density_field =
536 bc->block->template getData<DensityField>(m_density_field_id);
537
538 return {double_c(density_field->get(bc->cell))};
539 }
540
541 [[nodiscard]] std::vector<double>
543 Utils::Vector3i const &upper_corner) const override {
544 std::vector<double> out;
545 if (auto const ci = get_interval(lower_corner, upper_corner)) {
546 auto const &lattice = get_lattice();
547 auto const &block = *(lattice.get_blocks()->begin());
548 auto const density_field =
549 block.template getData<DensityField>(m_density_field_id);
550 auto const lower_cell = ci->min();
551 auto const upper_cell = ci->max();
552 auto const n_values = ci->numCells();
553 out.reserve(n_values);
554 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
555 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
556 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
557 out.emplace_back(density_field->get(Cell{x, y, z}));
558 }
559 }
560 }
561 assert(out.size() == n_values);
562 }
563 return out;
564 }
565
566 void set_slice_density(Utils::Vector3i const &lower_corner,
567 Utils::Vector3i const &upper_corner,
568 std::vector<double> const &density) override {
569 if (auto const ci = get_interval(lower_corner, upper_corner)) {
570 auto const &lattice = get_lattice();
571 auto &block = *(lattice.get_blocks()->begin());
572 auto density_field =
573 block.template getData<DensityField>(m_density_field_id);
574 auto it = density.begin();
575 auto const lower_cell = ci->min();
576 auto const upper_cell = ci->max();
577 assert(density.size() == ci->numCells());
578 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
579 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
580 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
581 density_field->get(Cell{x, y, z}) = FloatType_c(*it);
582 ++it;
583 }
584 }
585 }
586 }
587 }
588
589 void clear_flux_boundaries() override {
591 }
592
593 void clear_density_boundaries() override {
595 }
596
598 Utils::Vector3d const &flux) override {
599 auto bc = get_block_and_cell(get_lattice(), node, true);
600 if (!bc)
601 return false;
602
603 m_boundary_flux->set_node_value_at_boundary(
604 node, to_vector3<FloatType>(flux), *bc);
605
606 return true;
607 }
608
609 [[nodiscard]] std::optional<Utils::Vector3d>
611 bool consider_ghosts = false) const override {
612 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
613 if (!bc or !m_boundary_flux->node_is_boundary(node))
614 return std::nullopt;
615
616 return {to_vector3d(m_boundary_flux->get_node_value_at_boundary(node))};
617 }
618
620 auto bc = get_block_and_cell(get_lattice(), node, true);
621 if (!bc)
622 return false;
623
624 m_boundary_flux->remove_node_from_boundary(node, *bc);
625
626 return true;
627 }
628
630 double density) override {
631 auto bc = get_block_and_cell(get_lattice(), node, true);
632 if (!bc)
633 return false;
634
635 m_boundary_density->set_node_value_at_boundary(node, FloatType_c(density),
636 *bc);
637
638 return true;
639 }
640
641 [[nodiscard]] std::optional<double>
643 bool consider_ghosts = false) const override {
644 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
645 if (!bc or !m_boundary_density->node_is_boundary(node))
646 return std::nullopt;
647
648 return {double_c(m_boundary_density->get_node_value_at_boundary(node))};
649 }
650
652 Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
653 std::vector<std::optional<double>> const &density) override {
654 if (auto const ci = get_interval(lower_corner, upper_corner)) {
655 auto const &lattice = get_lattice();
656 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
657 auto const lower_cell = ci->min();
658 auto const upper_cell = ci->max();
659 auto it = density.begin();
660 assert(density.size() == ci->numCells());
661 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
662 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
663 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
664 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
665 auto const bc = get_block_and_cell(lattice, node, false);
666 auto const &opt = *it;
667 if (opt) {
668 m_boundary_density->set_node_value_at_boundary(
669 node, FloatType_c(*opt), *bc);
670 } else {
671 m_boundary_density->remove_node_from_boundary(node, *bc);
672 }
673 ++it;
674 }
675 }
676 }
677 }
678 }
679
680 [[nodiscard]] std::vector<std::optional<double>>
682 Utils::Vector3i const &lower_corner,
683 Utils::Vector3i const &upper_corner) const override {
684 std::vector<std::optional<double>> out;
685 if (auto const ci = get_interval(lower_corner, upper_corner)) {
686 auto const &lattice = get_lattice();
687 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
688 auto const lower_cell = ci->min();
689 auto const upper_cell = ci->max();
690 auto const n_values = ci->numCells();
691 out.reserve(n_values);
692 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
693 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
694 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
695 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
696 if (m_boundary_density->node_is_boundary(node)) {
697 out.emplace_back(double_c(
698 m_boundary_density->get_node_value_at_boundary(node)));
699 } else {
700 out.emplace_back(std::nullopt);
701 }
702 }
703 }
704 }
705 assert(out.size() == n_values);
706 }
707 return out;
708 }
709
711 Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
712 std::vector<std::optional<Utils::Vector3d>> const &flux) override {
713 if (auto const ci = get_interval(lower_corner, upper_corner)) {
714 auto const &lattice = get_lattice();
715 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
716 auto const lower_cell = ci->min();
717 auto const upper_cell = ci->max();
718 auto it = flux.begin();
719 assert(flux.size() == ci->numCells());
720 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
721 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
722 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
723 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
724 auto const bc = get_block_and_cell(lattice, node, false);
725 auto const &opt = *it;
726 if (opt) {
727 m_boundary_flux->set_node_value_at_boundary(
728 node, to_vector3<FloatType>(*opt), *bc);
729 } else {
730 m_boundary_flux->remove_node_from_boundary(node, *bc);
731 }
732 ++it;
733 }
734 }
735 }
736 }
737 }
738
739 [[nodiscard]] std::vector<std::optional<Utils::Vector3d>>
741 Utils::Vector3i const &lower_corner,
742 Utils::Vector3i const &upper_corner) const override {
743 std::vector<std::optional<Utils::Vector3d>> out;
744 if (auto const ci = get_interval(lower_corner, upper_corner)) {
745 auto const &lattice = get_lattice();
746 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
747 auto const lower_cell = ci->min();
748 auto const upper_cell = ci->max();
749 auto const n_values = ci->numCells();
750 out.reserve(n_values);
751 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
752 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
753 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
754 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
755 if (m_boundary_flux->node_is_boundary(node)) {
756 out.emplace_back(to_vector3d(
757 m_boundary_flux->get_node_value_at_boundary(node)));
758 } else {
759 out.emplace_back(std::nullopt);
760 }
761 }
762 }
763 }
764 assert(out.size() == n_values);
765 }
766 return out;
767 }
768
769 [[nodiscard]] std::vector<bool>
771 Utils::Vector3i const &upper_corner) const override {
772 std::vector<bool> out;
773 if (auto const ci = get_interval(lower_corner, upper_corner)) {
774 auto const &lattice = get_lattice();
775 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
776 auto const lower_cell = ci->min();
777 auto const upper_cell = ci->max();
778 auto const n_values = ci->numCells();
779 out.reserve(n_values);
780 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
781 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
782 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
783 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
784 out.emplace_back(m_boundary_density->node_is_boundary(node) or
785 m_boundary_flux->node_is_boundary(node));
786 }
787 }
788 }
789 assert(out.size() == n_values);
790 }
791 return out;
792 }
793
795 auto bc = get_block_and_cell(get_lattice(), node, true);
796 if (!bc)
797 return false;
798
799 m_boundary_density->remove_node_from_boundary(node, *bc);
800
801 return true;
802 }
803
804 [[nodiscard]] std::optional<bool>
806 bool consider_ghosts) const override {
807 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
808 if (!bc)
809 return std::nullopt;
810
811 return {m_boundary_flux->node_is_boundary(node)};
812 }
813
814 [[nodiscard]] std::optional<bool>
816 bool consider_ghosts) const override {
817 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
818 if (!bc)
819 return std::nullopt;
820
821 return {m_boundary_density->node_is_boundary(node)};
822 }
823
824 [[nodiscard]] std::optional<bool>
826 bool consider_ghosts = false) const override {
827 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
828 if (!bc)
829 return std::nullopt;
830
831 return {m_boundary_density->node_is_boundary(node) or
832 m_boundary_flux->node_is_boundary(node)};
833 }
834
836 const std::vector<int> &raster_flat,
837 const std::vector<double> &data_flat) override {
838 auto const grid_size = get_lattice().get_grid_dimensions();
839 auto const data = fill_3D_vector_array(data_flat, grid_size);
840 set_boundary_from_grid(*m_boundary_flux, get_lattice(), raster_flat, data);
842 }
843
845 const std::vector<int> &raster_flat,
846 const std::vector<double> &data_flat) override {
847 auto const grid_size = get_lattice().get_grid_dimensions();
848 auto const data = fill_3D_scalar_array(data_flat, grid_size);
850 data);
852 }
853
854 void reallocate_flux_boundary_field() { m_boundary_flux->boundary_update(); }
855
857 m_boundary_density->boundary_update();
858 }
859
860 [[nodiscard]] LatticeWalberla const &get_lattice() const noexcept override {
861 return *m_lattice;
862 }
863
864 void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override {
865 field::FlagFieldCellFilter<FlagField> fluid_filter(m_flag_field_density_id);
866 fluid_filter.addFlag(Boundary_flag);
867 vtk_obj.addCellExclusionFilter(fluid_filter);
868 }
869
870protected:
871 template <typename Field_T, uint_t F_SIZE_ARG, typename OutputType>
872 class VTKWriter : public vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG> {
873 public:
874 VTKWriter(ConstBlockDataID const &block_id, std::string const &id,
875 FloatType unit_conversion)
876 : vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG>(id),
877 m_block_id(block_id), m_field(nullptr),
878 m_conversion(unit_conversion) {}
879
880 protected:
881 void configure() override {
882 WALBERLA_ASSERT_NOT_NULLPTR(this->block_);
883 m_field = this->block_->template getData<Field_T>(m_block_id);
884 }
885
886 ConstBlockDataID const m_block_id;
887 Field_T const *m_field;
888 FloatType const m_conversion;
889 };
890
891 template <typename OutputType = float,
893 class DensityVTKWriter : public VTKWriter<DensityField, 1u, OutputType> {
894 public:
895 using VTKWriter<DensityField, 1u, OutputType>::VTKWriter;
896
897 protected:
898 OutputType evaluate(cell_idx_t const x, cell_idx_t const y,
899 cell_idx_t const z, cell_idx_t const) override {
900 WALBERLA_ASSERT_NOT_NULLPTR(this->m_field);
901 auto const density = VectorTrait<typename DensityField::value_type>::get(
902 this->m_field->get(x, y, z, 0), uint_c(0));
903 return numeric_cast<OutputType>(this->m_conversion * density);
904 }
905 };
906
907public:
908 void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj,
909 LatticeModel::units_map const &units,
910 int flag_observables) override {
911 if (flag_observables & static_cast<int>(EKOutputVTK::density)) {
912 auto const unit_conversion = FloatType_c(units.at("density"));
913 vtk_obj.addCellDataWriter(make_shared<DensityVTKWriter<float>>(
914 m_density_field_id, "density", unit_conversion));
915 }
916 }
917
918 ~EKinWalberlaImpl() override = default;
919};
920
921} // namespace walberla
Vector implementation and trait types for boost qvm interoperability.
Definition Cell.hpp:96
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.
std::pair< Utils::Vector3i, Utils::Vector3i > get_local_grid_range() const
walberla::blockforest::StructuredBlockForest Lattice_T
static DEVICE_QUALIFIER constexpr Vector< T, N > broadcast(typename Base::value_type const &value) noexcept
Create a vector that has all entries set to the same value.
Definition Vector.hpp:111
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
VTKWriter(ConstBlockDataID const &block_id, std::string const &id, FloatType unit_conversion)
Class that runs and controls the EK on waLBerla.
void update_density_boundary_from_shape(const std::vector< int > &raster_flat, const std::vector< double > &data_flat) override
double get_kT() const noexcept override
std::optional< CellInterval > get_interval(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const
blockforest::communication::UniformBufferedScheme< typename stencil::D3Q27 > FullCommunicator
std::shared_ptr< LatticeWalberla > m_lattice
Block forest.
bool get_friction_coupling() const noexcept override
void set_friction_coupling(bool friction_coupling) override
void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj, LatticeModel::units_map const &units, int flag_observables) override
std::vector< std::optional< double > > get_slice_density_at_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
void reset_density_boundary_handling(std::shared_ptr< BlockStorage > const &blocks)
void set_slice_density_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< std::optional< double > > const &density) override
~EKinWalberlaImpl() override=default
double get_valency() const noexcept override
std::optional< uint64_t > get_rng_state() const override
void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override
void set_rng_state(uint64_t counter) override
std::unique_ptr< DiffusiveFluxKernelElectrostatic > m_diffusive_flux_electrostatic
walberla::FlagField< walberla::uint8_t > FlagField
void set_slice_flux_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< std::optional< Utils::Vector3d > > const &flux) override
bool remove_node_from_flux_boundary(Utils::Vector3i const &node) override
bool remove_node_from_density_boundary(Utils::Vector3i const &node) override
std::optional< Utils::Vector3d > get_node_flux_at_boundary(Utils::Vector3i const &node, bool consider_ghosts=false) const override
void set_diffusion(double diffusion) override
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::size_t stencil_size() const noexcept override
void update_flux_boundary_from_shape(const std::vector< int > &raster_flat, const std::vector< double > &data_flat) override
bool set_node_density_boundary(Utils::Vector3i const &node, double density) override
bool is_thermalized() const noexcept override
void set_ext_efield(Utils::Vector3d const &field) override
std::vector< double > get_slice_density(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
std::unique_ptr< BoundaryModelDensity > m_boundary_density
std::optional< bool > get_node_is_flux_boundary(Utils::Vector3i const &node, bool consider_ghosts) const override
unsigned int get_seed() const noexcept override
bool get_advection() const noexcept override
std::size_t get_density_id() const noexcept override
GhostLayerField< FloatType, FluxCount > FluxField
double get_diffusion() 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
std::optional< double > get_node_density(Utils::Vector3i const &node, bool consider_ghosts=false) const override
std::shared_ptr< FullCommunicator > m_full_communication
std::unique_ptr< BoundaryModelFlux > m_boundary_flux
void set_slice_density(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner, std::vector< double > const &density) override
FlagUID const Domain_flag
Flag for domain cells, i.e.
bool is_double_precision() const noexcept override
void set_kT(double kT) override
void set_valency(double valency) override
std::optional< double > get_node_density_at_boundary(Utils::Vector3i const &node, bool consider_ghosts=false) const override
void integrate(std::size_t potential_id, std::size_t velocity_id, std::size_t force_id) override
LatticeWalberla::Lattice_T BlockStorage
std::optional< bool > get_node_is_density_boundary(Utils::Vector3i const &node, bool consider_ghosts) const override
GhostLayerField< FloatType, 1 > DensityField
Utils::Vector3d get_ext_efield() const noexcept override
LatticeWalberla const & get_lattice() const noexcept override
bool set_node_density(Utils::Vector3i const &node, double density) override
std::unique_ptr< ContinuityKernel > m_continuity
void clear_density_boundaries() override
std::optional< bool > get_node_is_boundary(Utils::Vector3i const &node, bool consider_ghosts=false) const override
std::unique_ptr< DiffusiveFluxKernel > m_diffusive_flux
void set_advection(bool advection) override
std::vector< std::optional< Utils::Vector3d > > get_slice_flux_at_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
void reset_flux_boundary_handling(std::shared_ptr< BlockStorage > const &blocks)
std::vector< bool > get_slice_is_boundary(Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner) const override
static double * block(double *p, std::size_t index, std::size_t size)
Definition elc.cpp:172
static FUNC_PREFIX double *RESTRICT 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 int64_t const uint32_t uint32_t uint32_t double double double double double uint32_t seed
\file PackInfoPdfDoublePrecision.cpp \author pystencils
std::optional< BlockAndCell > get_block_and_cell(::LatticeWalberla const &lattice, Utils::Vector3i const &node, bool consider_ghost_layers)
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::vector< Utils::Vector3d > fill_3D_vector_array(std::vector< double > const &vec_flat, Utils::Vector3i const &grid_size)
Definition boundary.hpp:36