ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
script_interface/walberla/VTKHandle.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2021-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 "config/config.hpp"
23
24#ifdef WALBERLA
25
28
31
32#include <boost/algorithm/string/join.hpp>
33
34#include <algorithm>
35#include <filesystem>
36#include <memory>
37#include <sstream>
38#include <stdexcept>
39#include <string>
40#include <unordered_map>
41#include <vector>
42
44
45template <class Field>
46class VTKHandleBase : public AutoParameters<VTKHandleBase<Field>> {
47private:
48 int m_delta_N;
49 int m_flag_obs;
50 std::string m_identifier;
51 std::string m_base_folder;
52 std::string m_prefix;
53 std::shared_ptr<::VTKHandle> m_vtk_handle;
54 std::weak_ptr<Field> m_field;
56 std::vector<Variant> m_pending_arguments;
57
58 [[nodiscard]] auto get_vtk_uid() const {
59 return m_base_folder + '/' + m_identifier;
60 }
61
62 [[nodiscard]] std::shared_ptr<Field> get_field_instance() const {
63 if (auto const field = m_field.lock()) {
64 return field;
65 }
66 auto const has_lattice_expired = m_pending_arguments.empty();
67 auto const err_expired = "Attempted access to an expired lattice object";
68 auto const err_detached = "This VTK object isn't attached to a lattice";
69 throw std::runtime_error(has_lattice_expired ? err_expired : err_detached);
70 }
71
72 virtual std::unordered_map<std::string, int> const &get_obs_map() const = 0;
73
74 [[nodiscard]] auto get_valid_observable_names() const {
75 std::vector<std::string> names{};
76 for (auto const &kv : get_obs_map()) {
77 names.emplace_back(kv.first);
78 }
79 std::sort(names.begin(), names.end());
80 return names;
81 }
82
83 [[nodiscard]] int
84 deserialize_obs_flag(std::vector<std::string> const &names) const {
85 int flag{0};
86 auto const &obs_map = get_obs_map();
87 for (auto const &name : names) {
88 if (obs_map.count(name) == 0) {
89 auto const valid_names = get_valid_observable_names();
90 std::stringstream message;
91 message << "Only the following VTK observables are supported: ['"
92 << boost::algorithm::join(valid_names, "', '") << "'], got '"
93 << name << "'";
94 throw std::invalid_argument(message.str());
95 }
96 flag |= obs_map.at(name);
97 }
98 return flag;
99 }
100
101 [[nodiscard]] Variant serialize_obs_flag(int flag) const {
102 std::vector<Variant> observables{};
103 for (auto const &kv : get_obs_map()) {
104 if (flag & kv.second) {
105 observables.emplace_back(kv.first);
106 }
107 }
108 return observables;
109 }
110
111public:
113 constexpr auto read_only = AutoParameter::read_only;
115 {"enabled", read_only, [this]() { return m_vtk_handle->enabled; }},
116 {"delta_N", read_only, [this]() { return m_delta_N; }},
117 {"vtk_uid", read_only, [this]() { return get_vtk_uid(); }},
118 {"identifier", read_only, [this]() { return m_identifier; }},
119 {"base_folder", read_only, [this]() { return m_base_folder; }},
120 {"prefix", read_only, [this]() { return m_prefix; }},
121 {"observables", read_only,
122 [this]() { return serialize_obs_flag(m_flag_obs); }},
123 {"execution_count", read_only,
124 [this]() { return m_vtk_handle->execution_count; }},
125 {"units", read_only,
126 [this]() { return make_unordered_map_of_variants(m_units); }},
127 });
128 }
129
130private:
131 void do_construct(VariantMap const &params) override {
132 m_delta_N = get_value<int>(params, "delta_N");
133 m_identifier = get_value<std::string>(params, "identifier");
134 m_base_folder = get_value<std::string>(params, "base_folder");
135 m_prefix = get_value<std::string>(params, "prefix");
136 auto const is_enabled = get_value<bool>(params, "enabled");
137 auto const execution_count = get_value<int>(params, "execution_count");
139 m_flag_obs = deserialize_obs_flag(
140 get_value<std::vector<std::string>>(params, "observables"));
141 if (m_delta_N < 0) {
142 throw std::domain_error("Parameter 'delta_N' must be >= 0");
143 }
144 if (m_identifier.empty()) {
145 throw std::domain_error("Parameter 'identifier' cannot be empty");
146 }
147 if (m_identifier.find(std::filesystem::path::preferred_separator) !=
148 std::string::npos) {
149 throw std::invalid_argument(
150 "Parameter 'identifier' cannot be a filepath");
151 }
152 });
153 m_pending_arguments.emplace_back(is_enabled);
154 m_pending_arguments.emplace_back(execution_count);
155 }
156
157protected:
158 Variant do_call_method(std::string const &name,
159 VariantMap const &params) override {
160 if (name == "enable") {
162 if (m_delta_N == 0) {
163 throw std::runtime_error("Manual VTK callbacks cannot be enabled");
164 }
165 get_field_instance()->switch_vtk(get_vtk_uid(), true);
166 });
167 return {};
168 }
169 if (name == "disable") {
171 if (m_delta_N == 0) {
172 throw std::runtime_error("Manual VTK callbacks cannot be disabled");
173 }
174 get_field_instance()->switch_vtk(get_vtk_uid(), false);
175 });
176 return {};
177 }
178 if (name == "write") {
180 if (m_delta_N) {
181 throw std::runtime_error("Automatic VTK callbacks cannot be "
182 "triggered manually");
183 }
184 get_field_instance()->write_vtk(get_vtk_uid());
185 });
186 return {};
187 }
188 if (name == "get_valid_observable_names") {
189 return make_vector_of_variants(get_valid_observable_names());
190 }
191
192 return {};
193 }
194
195public:
196 void detach_from_lattice() { m_field.reset(); }
197
198 void attach_to_lattice(std::weak_ptr<Field> field,
199 ::LatticeModel::units_map const &units) {
200 auto const was_attached_once = m_pending_arguments.empty();
201 if (not m_field.expired()) {
202 throw std::runtime_error("Cannot attach VTK object to multiple lattices");
203 }
204 if (was_attached_once) {
205 throw std::runtime_error("Detached VTK objects cannot be attached again");
206 }
207 assert(m_pending_arguments.size() == 2u);
208 auto const is_enabled = get_value<bool>(m_pending_arguments[0]);
209 auto const execution_count = get_value<int>(m_pending_arguments[1]);
210 m_units = units;
211 m_field = field;
212 auto instance = get_field_instance();
213 m_vtk_handle =
214 instance->create_vtk(m_delta_N, execution_count, m_flag_obs, m_units,
215 m_identifier, m_base_folder, m_prefix);
216 if (m_delta_N and not is_enabled) {
217 instance->switch_vtk(get_vtk_uid(), false);
218 }
219 m_pending_arguments.clear();
220 }
221};
222
223} // namespace ScriptInterface::walberla
224
225#endif // WALBERLA
std::unordered_map< std::string, double > units_map
Bind parameters in the script interface.
void add_parameters(std::vector< AutoParameter > &&params)
virtual void parallel_try_catch(std::function< void()> const &cb) const =0
boost::string_ref name() const
Context * context() const
Responsible context.
void attach_to_lattice(std::weak_ptr< Field > field, ::LatticeModel::units_map const &units)
virtual std::unordered_map< std::string, int > const & get_obs_map() const =0
Variant do_call_method(std::string const &name, VariantMap const &params) override
This file contains the defaults for ESPResSo.
T get_value(Variant const &v)
Extract value of specific type T from a Variant.
std::unordered_map< std::string, Variant > VariantMap
Definition Variant.hpp:69
auto make_unordered_map_of_variants(std::unordered_map< K, V > const &v)
Definition Variant.hpp:80
auto make_vector_of_variants(std::vector< T > const &v)
Definition Variant.hpp:88
boost::make_recursive_variant< None, bool, int, std::size_t, double, std::string, ObjectRef, Utils::Vector3b, Utils::Vector3i, Utils::Vector2d, Utils::Vector3d, Utils::Vector4d, std::vector< int >, std::vector< double >, std::vector< boost::recursive_variant_ >, std::unordered_map< int, boost::recursive_variant_ >, std::unordered_map< std::string, boost::recursive_variant_ > >::type Variant
Possible types for parameters.
Definition Variant.hpp:67
static SteepestDescentParameters params
Currently active steepest descent instance.
static constexpr const ReadOnly read_only