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
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 <vector>
51
52namespace walberla {
53
54/** @brief Class that runs and controls the EK on waLBerla. */
55template <std::size_t FluxCount = 13, typename FloatType = double>
57 using ContinuityKernel =
59 using DiffusiveFluxKernel =
61 using AdvectiveFluxKernel =
63 using FrictionCouplingKernel =
65 using DiffusiveFluxKernelElectrostatic =
67
68 using Dirichlet = typename detail::KernelTrait<FloatType>::Dirichlet;
69 using FixedFlux = typename detail::KernelTrait<FloatType>::FixedFlux;
70
71protected:
72 // Type definitions
73 using FluxField = GhostLayerField<FloatType, FluxCount>;
74 using FlagField = walberla::FlagField<walberla::uint8_t>;
75 using DensityField = GhostLayerField<FloatType, 1>;
76
79
80public:
81 template <typename T> FloatType FloatType_c(T t) {
82 return numeric_cast<FloatType>(t);
83 }
84
85 [[nodiscard]] std::size_t stencil_size() const noexcept override {
86 return FluxCount;
87 }
88
89 [[nodiscard]] bool is_double_precision() const noexcept override {
90 return std::is_same<FloatType, double>::value;
91 }
92
93private:
94 FloatType m_diffusion;
95 FloatType m_kT;
96 FloatType m_valency;
97 Utils::Vector3d m_ext_efield;
98 bool m_advection;
99 bool m_friction_coupling;
100
101protected:
102 // Block data access handles
105
106 BlockDataID m_flux_field_id;
108
111
112 /** Flag for domain cells, i.e. all cells. */
113 FlagUID const Domain_flag{"domain"};
114 /** Flag for boundary cells. */
115 FlagUID const Boundary_flag{"boundary"};
116
117 /** Block forest */
118 std::shared_ptr<LatticeWalberla> m_lattice;
119
120 std::unique_ptr<BoundaryModelDensity> m_boundary_density;
121 std::unique_ptr<BoundaryModelFlux> m_boundary_flux;
122
123 std::unique_ptr<ContinuityKernel> m_continuity;
124
125 // ResetFlux + external force
126 // TODO: kernel for that
127 // std::shared_ptr<ResetForce<PdfField, VectorField>> m_reset_force;
128
129 [[nodiscard]] std::optional<CellInterval>
130 get_interval(Utils::Vector3i const &lower_corner,
131 Utils::Vector3i const &upper_corner) const {
132 auto const &lattice = get_lattice();
133 auto const &cell_min = lower_corner;
134 auto const cell_max = upper_corner - Utils::Vector3i::broadcast(1);
135 auto const lower_bc = get_block_and_cell(lattice, cell_min, true);
136 auto const upper_bc = get_block_and_cell(lattice, cell_max, true);
137 if (not lower_bc or not upper_bc) {
138 return std::nullopt;
139 }
140 assert(&(*(lower_bc->block)) == &(*(upper_bc->block)));
141 return {CellInterval(lower_bc->cell, upper_bc->cell)};
142 }
143
145 auto const &blocks = get_lattice().get_blocks();
146 m_boundary_density = std::make_unique<BoundaryModelDensity>(
148 }
149
151 auto const &blocks = get_lattice().get_blocks();
152 m_boundary_flux = std::make_unique<BoundaryModelFlux>(
154 }
155
156 using FullCommunicator = blockforest::communication::UniformBufferedScheme<
157 typename stencil::D3Q27>;
158 std::shared_ptr<FullCommunicator> m_full_communication;
159
160public:
161 EKinWalberlaImpl(std::shared_ptr<LatticeWalberla> lattice, double diffusion,
162 double kT, double valency, Utils::Vector3d const &ext_efield,
163 double density, bool advection, bool friction_coupling)
164 : m_diffusion(FloatType_c(diffusion)), m_kT(FloatType_c(kT)),
165 m_valency(FloatType_c(valency)), m_ext_efield(ext_efield),
166 m_advection(advection), m_friction_coupling(friction_coupling),
167 m_lattice(std::move(lattice)) {
168 m_density_field_id = field::addToStorage<DensityField>(
169 m_lattice->get_blocks(), "density field", FloatType_c(density),
170 field::fzyx, m_lattice->get_ghost_layers());
172 field::addFlattenedShallowCopyToStorage<DensityField>(
173 m_lattice->get_blocks(), m_density_field_id,
174 "flattened density field");
175 m_flux_field_id = field::addToStorage<FluxField>(
176 m_lattice->get_blocks(), "flux field", FloatType{0}, field::fzyx,
177 m_lattice->get_ghost_layers());
179 field::addFlattenedShallowCopyToStorage<FluxField>(
180 m_lattice->get_blocks(), m_flux_field_id, "flattened flux field");
181
182 m_continuity = std::make_unique<ContinuityKernel>(
184
185 // Init boundary related stuff
186 m_flag_field_density_id = field::addFlagFieldToStorage<FlagField>(
187 m_lattice->get_blocks(), "flag field density",
188 m_lattice->get_ghost_layers());
190
191 m_flag_field_flux_id = field::addFlagFieldToStorage<FlagField>(
192 m_lattice->get_blocks(), "flag field flux",
193 m_lattice->get_ghost_layers());
195
197 std::make_shared<FullCommunicator>(m_lattice->get_blocks());
198 m_full_communication->addPackInfo(
199 std::make_shared<field::communication::PackInfo<DensityField>>(
201
202 // Synchronize ghost layers
204 }
205
206 // Global parameters
207 [[nodiscard]] double get_diffusion() const noexcept override {
208 return m_diffusion;
209 }
210 [[nodiscard]] double get_kT() const noexcept override { return m_kT; }
211 [[nodiscard]] double get_valency() const noexcept override {
212 return m_valency;
213 }
214 [[nodiscard]] bool get_advection() const noexcept override {
215 return m_advection;
216 }
217 [[nodiscard]] bool get_friction_coupling() const noexcept override {
218 return m_friction_coupling;
219 }
220 [[nodiscard]] Utils::Vector3d get_ext_efield() const noexcept override {
221 return m_ext_efield;
222 }
223
224 void set_diffusion(double diffusion) override {
225 m_diffusion = FloatType_c(diffusion);
226 }
227 void set_kT(double kT) override { m_kT = FloatType_c(kT); }
228 void set_valency(double valency) override {
229 m_valency = FloatType_c(valency);
230 }
231 void set_advection(bool advection) override { m_advection = advection; }
232 void set_friction_coupling(bool friction_coupling) override {
233 m_friction_coupling = friction_coupling;
234 }
235 void set_ext_efield(Utils::Vector3d const &field) override {
236 m_ext_efield = field;
237 }
238
239 void ghost_communication() override { (*m_full_communication)(); }
240
241private:
242 void kernel_boundary_density() {
243 for (auto &block : *m_lattice->get_blocks()) {
244 (*m_boundary_density)(&block);
245 }
246 }
247
248 void kernel_boundary_flux() {
249 for (auto &block : *m_lattice->get_blocks()) {
250 (*m_boundary_flux)(&block);
251 }
252 }
253
254 void kernel_continuity() {
255 for (auto &block : *m_lattice->get_blocks()) {
256 (*m_continuity).run(&block);
257 }
258 }
259
260 void kernel_diffusion() {
261 auto kernel = DiffusiveFluxKernel(m_flux_field_flattened_id,
264
265 for (auto &block : *m_lattice->get_blocks()) {
266 kernel.run(&block);
267 }
268 }
269
270 void kernel_advection(const std::size_t &velocity_id) {
271 auto kernel =
273 BlockDataID(velocity_id));
274 for (auto &block : *m_lattice->get_blocks()) {
275 kernel.run(&block);
276 }
277 }
278
279 void kernel_friction_coupling(const std::size_t &force_id) {
280 auto kernel = FrictionCouplingKernel(
281 BlockDataID(force_id), m_flux_field_flattened_id,
283 for (auto &block : *m_lattice->get_blocks()) {
284 kernel.run(&block);
285 }
286 }
287
288 void kernel_diffusion_electrostatic(const std::size_t &potential_id) {
289 auto const ext_field = get_ext_efield();
290 auto kernel = DiffusiveFluxKernelElectrostatic(
291 m_flux_field_flattened_id, BlockDataID(potential_id),
293 FloatType_c(ext_field[0]), FloatType_c(ext_field[1]),
294 FloatType_c(ext_field[2]), FloatType_c(get_kT()),
296 for (auto &block : *m_lattice->get_blocks()) {
297 kernel.run(&block);
298 }
299 }
300
301 void kernel_migration() {}
302
303 void updated_boundary_fields() {
304 m_boundary_flux->boundary_update();
305 m_boundary_density->boundary_update();
306 }
307
308protected:
309 void integrate_vtk_writers() override {
310 for (auto const &it : m_vtk_auto) {
311 auto &vtk_handle = it.second;
312 if (vtk_handle->enabled) {
313 vtk::writeFiles(vtk_handle->ptr)();
314 vtk_handle->execution_count++;
315 }
316 }
317 }
318
319public:
320 void integrate(std::size_t potential_id, std::size_t velocity_id,
321 std::size_t force_id) override {
322
323 updated_boundary_fields();
324
325 if (get_diffusion() == 0.)
326 return;
327
328 if (get_valency() != 0.) {
329 if (potential_id == walberla::BlockDataID{}) {
330 throw std::runtime_error("Walberla EK: electrostatic potential enabled "
331 "but no field accessible. potential id is " +
332 std::to_string(potential_id));
333 }
334 kernel_diffusion_electrostatic(potential_id);
335 } else {
336 kernel_diffusion();
337 }
338
339 kernel_migration();
340 kernel_boundary_flux();
341 // friction coupling
342 if (get_friction_coupling()) {
343 if (force_id == walberla::BlockDataID{}) {
344 throw std::runtime_error("Walberla EK: friction coupling enabled but "
345 "no force field accessible. force_id is " +
346 std::to_string(force_id) +
347 ". Hint: LB may be inactive.");
348 }
349 kernel_friction_coupling(force_id);
350 }
351
352 if (get_advection()) {
353 if (velocity_id == walberla::BlockDataID{}) {
354 throw std::runtime_error("Walberla EK: advection enabled but no "
355 "velocity field accessible. velocity_id is " +
356 std::to_string(velocity_id) +
357 ". Hint: LB may be inactive.");
358 }
359 kernel_advection(velocity_id);
360 }
361 kernel_continuity();
362
363 // is this the expected behavior when reactions are included?
364 kernel_boundary_density();
365
366 // Handle VTK writers
368 }
369
370 [[nodiscard]] std::size_t get_density_id() const noexcept override {
371 static_assert(std::is_same_v<std::size_t, walberla::uint_t>);
372 return static_cast<std::size_t>(m_density_field_id);
373 }
374
375 bool set_node_density(Utils::Vector3i const &node, double density) override {
376 auto bc = get_block_and_cell(get_lattice(), node, false);
377 if (!bc)
378 return false;
379
380 auto density_field =
381 bc->block->template getData<DensityField>(m_density_field_id);
382 density_field->get(bc->cell) = FloatType_c(density);
383
384 return true;
385 }
386
387 [[nodiscard]] std::optional<double>
389 bool consider_ghosts = false) const override {
390 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
391
392 if (!bc)
393 return std::nullopt;
394
395 auto const density_field =
396 bc->block->template getData<DensityField>(m_density_field_id);
397
398 return {double_c(density_field->get(bc->cell))};
399 }
400
401 [[nodiscard]] std::vector<double>
403 Utils::Vector3i const &upper_corner) const override {
404 std::vector<double> out;
405 if (auto const ci = get_interval(lower_corner, upper_corner)) {
406 auto const &lattice = get_lattice();
407 auto const &block = *(lattice.get_blocks()->begin());
408 auto const density_field =
409 block.template getData<DensityField>(m_density_field_id);
410 auto const lower_cell = ci->min();
411 auto const upper_cell = ci->max();
412 auto const n_values = ci->numCells();
413 out.reserve(n_values);
414 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
415 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
416 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
417 out.emplace_back(density_field->get(Cell{x, y, z}));
418 }
419 }
420 }
421 assert(out.size() == n_values);
422 }
423 return out;
424 }
425
426 void set_slice_density(Utils::Vector3i const &lower_corner,
427 Utils::Vector3i const &upper_corner,
428 std::vector<double> const &density) override {
429 if (auto const ci = get_interval(lower_corner, upper_corner)) {
430 auto const &lattice = get_lattice();
431 auto &block = *(lattice.get_blocks()->begin());
432 auto density_field =
433 block.template getData<DensityField>(m_density_field_id);
434 auto it = density.begin();
435 auto const lower_cell = ci->min();
436 auto const upper_cell = ci->max();
437 assert(density.size() == ci->numCells());
438 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
439 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
440 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
441 density_field->get(Cell{x, y, z}) = FloatType_c(*it);
442 ++it;
443 }
444 }
445 }
446 }
447 }
448
450
454
456 Utils::Vector3d const &flux) override {
457 auto bc = get_block_and_cell(get_lattice(), node, true);
458 if (!bc)
459 return false;
460
461 m_boundary_flux->set_node_value_at_boundary(
462 node, to_vector3<FloatType>(flux), *bc);
463
464 return true;
465 }
466
467 [[nodiscard]] std::optional<Utils::Vector3d>
469 bool consider_ghosts = false) const override {
470 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
471 if (!bc or !m_boundary_flux->node_is_boundary(node))
472 return std::nullopt;
473
474 return {to_vector3d(m_boundary_flux->get_node_value_at_boundary(node))};
475 }
476
478 auto bc = get_block_and_cell(get_lattice(), node, true);
479 if (!bc)
480 return false;
481
482 m_boundary_flux->remove_node_from_boundary(node, *bc);
483
484 return true;
485 }
486
488 double density) override {
489 auto bc = get_block_and_cell(get_lattice(), node, true);
490 if (!bc)
491 return false;
492
493 m_boundary_density->set_node_value_at_boundary(node, FloatType_c(density),
494 *bc);
495
496 return true;
497 }
498
499 [[nodiscard]] std::optional<double>
501 bool consider_ghosts = false) const override {
502 auto const bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
503 if (!bc or !m_boundary_density->node_is_boundary(node))
504 return std::nullopt;
505
506 return {double_c(m_boundary_density->get_node_value_at_boundary(node))};
507 }
508
510 Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
511 std::vector<std::optional<double>> const &density) override {
512 if (auto const ci = get_interval(lower_corner, upper_corner)) {
513 auto const &lattice = get_lattice();
514 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
515 auto const lower_cell = ci->min();
516 auto const upper_cell = ci->max();
517 auto it = density.begin();
518 assert(density.size() == ci->numCells());
519 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
520 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
521 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
522 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
523 auto const bc = get_block_and_cell(lattice, node, false);
524 auto const &opt = *it;
525 if (opt) {
526 m_boundary_density->set_node_value_at_boundary(
527 node, FloatType_c(*opt), *bc);
528 } else {
529 m_boundary_density->remove_node_from_boundary(node, *bc);
530 }
531 ++it;
532 }
533 }
534 }
535 }
536 }
537
538 [[nodiscard]] std::vector<std::optional<double>>
540 Utils::Vector3i const &lower_corner,
541 Utils::Vector3i const &upper_corner) const override {
542 std::vector<std::optional<double>> out;
543 if (auto const ci = get_interval(lower_corner, upper_corner)) {
544 auto const &lattice = get_lattice();
545 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
546 auto const lower_cell = ci->min();
547 auto const upper_cell = ci->max();
548 auto const n_values = ci->numCells();
549 out.reserve(n_values);
550 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
551 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
552 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
553 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
554 if (m_boundary_density->node_is_boundary(node)) {
555 out.emplace_back(double_c(
556 m_boundary_density->get_node_value_at_boundary(node)));
557 } else {
558 out.emplace_back(std::nullopt);
559 }
560 }
561 }
562 }
563 assert(out.size() == n_values);
564 }
565 return out;
566 }
567
569 Utils::Vector3i const &lower_corner, Utils::Vector3i const &upper_corner,
570 std::vector<std::optional<Utils::Vector3d>> const &flux) override {
571 if (auto const ci = get_interval(lower_corner, upper_corner)) {
572 auto const &lattice = get_lattice();
573 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
574 auto const lower_cell = ci->min();
575 auto const upper_cell = ci->max();
576 auto it = flux.begin();
577 assert(flux.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 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
582 auto const bc = get_block_and_cell(lattice, node, false);
583 auto const &opt = *it;
584 if (opt) {
585 m_boundary_flux->set_node_value_at_boundary(
586 node, to_vector3<FloatType>(*opt), *bc);
587 } else {
588 m_boundary_flux->remove_node_from_boundary(node, *bc);
589 }
590 ++it;
591 }
592 }
593 }
594 }
595 }
596
597 [[nodiscard]] std::vector<std::optional<Utils::Vector3d>>
599 Utils::Vector3i const &lower_corner,
600 Utils::Vector3i const &upper_corner) const override {
601 std::vector<std::optional<Utils::Vector3d>> out;
602 if (auto const ci = get_interval(lower_corner, upper_corner)) {
603 auto const &lattice = get_lattice();
604 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
605 auto const lower_cell = ci->min();
606 auto const upper_cell = ci->max();
607 auto const n_values = ci->numCells();
608 out.reserve(n_values);
609 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
610 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
611 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
612 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
613 if (m_boundary_flux->node_is_boundary(node)) {
614 out.emplace_back(to_vector3d(
615 m_boundary_flux->get_node_value_at_boundary(node)));
616 } else {
617 out.emplace_back(std::nullopt);
618 }
619 }
620 }
621 }
622 assert(out.size() == n_values);
623 }
624 return out;
625 }
626
627 [[nodiscard]] std::vector<bool>
629 Utils::Vector3i const &upper_corner) const override {
630 std::vector<bool> out;
631 if (auto const ci = get_interval(lower_corner, upper_corner)) {
632 auto const &lattice = get_lattice();
633 auto const local_offset = std::get<0>(lattice.get_local_grid_range());
634 auto const lower_cell = ci->min();
635 auto const upper_cell = ci->max();
636 auto const n_values = ci->numCells();
637 out.reserve(n_values);
638 for (auto x = lower_cell.x(); x <= upper_cell.x(); ++x) {
639 for (auto y = lower_cell.y(); y <= upper_cell.y(); ++y) {
640 for (auto z = lower_cell.z(); z <= upper_cell.z(); ++z) {
641 auto const node = local_offset + Utils::Vector3i{{x, y, z}};
642 out.emplace_back(m_boundary_density->node_is_boundary(node) or
643 m_boundary_flux->node_is_boundary(node));
644 }
645 }
646 }
647 assert(out.size() == n_values);
648 }
649 return out;
650 }
651
653 auto bc = get_block_and_cell(get_lattice(), node, true);
654 if (!bc)
655 return false;
656
657 m_boundary_density->remove_node_from_boundary(node, *bc);
658
659 return true;
660 }
661
662 [[nodiscard]] std::optional<bool>
664 bool consider_ghosts) const override {
665 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
666 if (!bc)
667 return std::nullopt;
668
669 return {m_boundary_flux->node_is_boundary(node)};
670 }
671
672 [[nodiscard]] std::optional<bool>
674 bool consider_ghosts) const override {
675 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
676 if (!bc)
677 return std::nullopt;
678
679 return {m_boundary_density->node_is_boundary(node)};
680 }
681
682 [[nodiscard]] std::optional<bool>
684 bool consider_ghosts = false) const override {
685 auto bc = get_block_and_cell(get_lattice(), node, consider_ghosts);
686 if (!bc)
687 return std::nullopt;
688
689 return {m_boundary_density->node_is_boundary(node) or
690 m_boundary_flux->node_is_boundary(node)};
691 }
692
694 const std::vector<int> &raster_flat,
695 const std::vector<double> &data_flat) override {
696 auto const grid_size = get_lattice().get_grid_dimensions();
697 auto const data = fill_3D_vector_array(data_flat, grid_size);
698 set_boundary_from_grid(*m_boundary_flux, get_lattice(), raster_flat, data);
700 }
701
703 const std::vector<int> &raster_flat,
704 const std::vector<double> &data_flat) override {
705 auto const grid_size = get_lattice().get_grid_dimensions();
706 auto const data = fill_3D_scalar_array(data_flat, grid_size);
708 data);
710 }
711
712 void reallocate_flux_boundary_field() { m_boundary_flux->boundary_update(); }
713
715 m_boundary_density->boundary_update();
716 }
717
718 [[nodiscard]] LatticeWalberla const &get_lattice() const noexcept override {
719 return *m_lattice;
720 }
721
722 void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override {
723 field::FlagFieldCellFilter<FlagField> fluid_filter(m_flag_field_density_id);
724 fluid_filter.addFlag(Boundary_flag);
725 vtk_obj.addCellExclusionFilter(fluid_filter);
726 }
727
728protected:
729 template <typename Field_T, uint_t F_SIZE_ARG, typename OutputType>
730 class VTKWriter : public vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG> {
731 public:
732 VTKWriter(ConstBlockDataID const &block_id, std::string const &id,
733 FloatType unit_conversion)
734 : vtk::BlockCellDataWriter<OutputType, F_SIZE_ARG>(id),
735 m_block_id(block_id), m_field(nullptr),
736 m_conversion(unit_conversion) {}
737
738 protected:
739 void configure() override {
740 WALBERLA_ASSERT_NOT_NULLPTR(this->block_);
741 m_field = this->block_->template getData<Field_T>(m_block_id);
742 }
743
744 ConstBlockDataID const m_block_id;
745 Field_T const *m_field;
746 FloatType const m_conversion;
747 };
748
749 template <typename OutputType = float,
751 class DensityVTKWriter : public VTKWriter<DensityField, 1u, OutputType> {
752 public:
753 using VTKWriter<DensityField, 1u, OutputType>::VTKWriter;
754
755 protected:
756 OutputType evaluate(cell_idx_t const x, cell_idx_t const y,
757 cell_idx_t const z, cell_idx_t const) override {
758 WALBERLA_ASSERT_NOT_NULLPTR(this->m_field);
759 auto const density = VectorTrait<typename DensityField::value_type>::get(
760 this->m_field->get(x, y, z, 0), uint_c(0));
761 return numeric_cast<OutputType>(this->m_conversion * density);
762 }
763 };
764
765public:
766 void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj,
767 LatticeModel::units_map const &units,
768 int flag_observables) override {
769 if (flag_observables & static_cast<int>(EKOutputVTK::density)) {
770 auto const unit_conversion = FloatType_c(units.at("density"));
771 vtk_obj.addCellDataWriter(make_shared<DensityVTKWriter<float>>(
772 m_density_field_id, "density", unit_conversion));
773 }
774 }
775
776 ~EKinWalberlaImpl() override = default;
777};
778
779} // namespace walberla
Vector implementation and trait types for boost qvm interoperability.
__shared__ int node[MAXDEPTH *THREADS5/WARPSIZE]
Definition Cell.hpp:98
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.
static Vector< T, N > broadcast(T const &s)
Create a vector that has all entries set to one value.
Definition Vector.hpp:110
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 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
void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) override
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)
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
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
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
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
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
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:174
std::optional< BlockAndCell > get_block_and_cell(::LatticeWalberla const &lattice, Utils::Vector3i const &node, bool consider_ghost_layers)
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
Utils::Vector3d to_vector3d(Vector3< float > const &v)