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 ESPRESSO_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_obs_flag;
50 bool m_force_pvtu;
51 std::string m_identifier;
52 std::filesystem::path m_base_folder;
53 std::string m_prefix;
54 std::shared_ptr<::VTKHandle> m_vtk_handle;
55 std::weak_ptr<Field> m_field;
57 std::vector<Variant> m_pending_arguments;
58
59 [[nodiscard]] auto get_vtk_uid() const {
60 return m_base_folder.generic_string() + "/" + m_identifier;
61 }
62
63 [[nodiscard]] std::shared_ptr<Field> get_field_instance() const {
64 if (auto const field = m_field.lock()) {
65 return field;
66 }
67 auto const has_lattice_expired = m_pending_arguments.empty();
68 auto const err_expired = "Attempted access to an expired lattice object";
69 auto const err_detached = "This VTK object isn't attached to a lattice";
70 throw std::runtime_error(has_lattice_expired ? err_expired : err_detached);
71 }
72
73 virtual std::unordered_map<std::string, int> const &get_obs_map() const = 0;
74
75 [[nodiscard]] auto get_valid_observable_names() const {
76 auto const view = std::views::elements<0>(get_obs_map());
77 std::vector<std::string> names{view.begin(), view.end()};
78 std::ranges::sort(names);
79 return names;
80 }
81
82 [[nodiscard]] int
83 deserialize_obs_flag(std::vector<std::string> const &names) const {
84 int obs_flag{0};
85 auto const &obs_map = get_obs_map();
86 for (auto const &name : names) {
87 if (not obs_map.contains(name)) {
88 auto const valid_names = get_valid_observable_names();
89 std::stringstream message;
90 message << "Only the following VTK observables are supported: ['"
91 << boost::algorithm::join(valid_names, "', '") << "'], got '"
92 << name << "'";
93 throw std::invalid_argument(message.str());
94 }
95 obs_flag |= obs_map.at(name);
96 }
97 return obs_flag;
98 }
99
100 [[nodiscard]] Variant serialize_obs_flag(int obs_flag) const {
101 std::vector<Variant> observables{};
102 for (auto const &[name, flag] : get_obs_map()) {
103 if (flag & obs_flag) {
104 observables.emplace_back(name);
105 }
106 }
107 return observables;
108 }
109
110public:
112 constexpr auto read_only = AutoParameter::read_only;
114 {"enabled", read_only, [this]() { return m_vtk_handle->enabled; }},
115 {"delta_N", read_only, [this]() { return m_delta_N; }},
116 {"vtk_uid", read_only, [this]() { return get_vtk_uid(); }},
117 {"identifier", read_only, [this]() { return m_identifier; }},
118 {"base_folder", read_only, [this]() { return m_base_folder; }},
119 {"prefix", read_only, [this]() { return m_prefix; }},
120 {"force_pvtu", read_only, [this]() { return m_force_pvtu; }},
121 {"observables", read_only,
122 [this]() { return serialize_obs_flag(m_obs_flag); }},
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::filesystem::path>(params, "base_folder");
135 m_prefix = get_value<std::string>(params, "prefix");
136 m_force_pvtu = get_value_or<bool>(params, "force_pvtu", false);
137 auto const is_enabled = get_value<bool>(params, "enabled");
138 auto const execution_count = get_value<int>(params, "execution_count");
140 m_obs_flag = deserialize_obs_flag(
141 get_value<std::vector<std::string>>(params, "observables"));
142 if (m_delta_N < 0) {
143 throw std::domain_error("Parameter 'delta_N' must be >= 0");
144 }
145 if (m_identifier.empty()) {
146 throw std::domain_error("Parameter 'identifier' cannot be empty");
147 }
148 if (m_identifier.find(std::filesystem::path::preferred_separator) !=
149 std::string::npos) {
150 throw std::invalid_argument(
151 "Parameter 'identifier' cannot be a filepath");
152 }
153 });
154 m_pending_arguments.emplace_back(is_enabled);
155 m_pending_arguments.emplace_back(execution_count);
156 }
157
158protected:
159 Variant do_call_method(std::string const &name, VariantMap const &) 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,
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 auto const base_folder = m_base_folder.generic_string();
211 m_units = units;
212 m_field = field;
213 auto instance = get_field_instance();
214 m_vtk_handle =
215 instance->create_vtk(m_delta_N, execution_count, m_obs_flag, m_units,
216 m_identifier, base_folder, m_prefix, m_force_pvtu);
217 if (m_delta_N and not is_enabled) {
218 instance->switch_vtk(get_vtk_uid(), false);
219 }
220 m_pending_arguments.clear();
221 }
222};
223
224} // namespace ScriptInterface::walberla
225
226#endif // ESPRESSO_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
Context * context() const
Responsible context.
std::string_view name() const
Variant do_call_method(std::string const &name, VariantMap const &) override
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
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:133
auto make_unordered_map_of_variants(std::unordered_map< K, V > const &v)
Definition Variant.hpp:144
auto make_vector_of_variants(std::vector< T > const &v)
Definition Variant.hpp:148
make_recursive_variant< ObjectRef > Variant
Possible types for parameters.
Definition Variant.hpp:131
static SteepestDescentParameters params
Currently active steepest descent instance.
static constexpr const ReadOnly read_only
Recursive variant implementation.
Definition Variant.hpp:84