ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
script_interface/walberla/LatticeModel.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 "LatticeWalberla.hpp"
23
25
29
30#include <algorithm>
31#include <memory>
32#include <stdexcept>
33#include <string>
34#include <vector>
35
37
38template <class Method, class VTKHandle>
40 : public AutoParameters<LatticeModel<Method, VTKHandle>, System::Leaf> {
41protected:
42 std::shared_ptr<LatticeWalberla> m_lattice;
43 std::shared_ptr<Method> m_instance;
44 std::vector<std::shared_ptr<VTKHandle>> m_vtk_writers;
45
46 virtual ::LatticeModel::units_map
48 virtual void make_instance(VariantMap const &params) = 0;
49
50 auto find_vtk(std::shared_ptr<VTKHandle> const &vtk) const {
51 return std::ranges::find(m_vtk_writers, vtk);
52 }
53
57
58public:
60 VariantMap const &params) override {
61 if (method_name == "add_vtk_writer") {
62 auto vtk = get_value<std::shared_ptr<VTKHandle>>(params, "vtk");
63 auto const needle = find_vtk(vtk);
65 if (needle != m_vtk_writers.end()) {
66 throw std::runtime_error(
67 "VTK object is already attached to this lattice");
68 }
69 vtk->attach_to_lattice(m_instance,
71 m_vtk_writers.emplace_back(vtk);
72 });
73 return {};
74 }
75 if (method_name == "remove_vtk_writer") {
76 auto const vtk = get_value<std::shared_ptr<VTKHandle>>(params, "vtk");
77 auto const needle = find_vtk(vtk);
79 if (needle == m_vtk_writers.end()) {
80 throw std::runtime_error(
81 "VTK object is not attached to this lattice");
82 }
83 vtk->detach_from_lattice();
84 });
85 m_vtk_writers.erase(needle);
86 return {};
87 }
88 if (method_name == "clear_vtk_writers") {
89 for (auto const &vtk : m_vtk_writers) {
90 vtk->detach_from_lattice();
91 }
92 m_vtk_writers.clear();
93 }
94 return {};
95 }
96};
97
98} // namespace ScriptInterface::walberla
Bind parameters in the script interface.
virtual void parallel_try_catch(std::function< void()> const &cb) const =0
Context * context() const
Responsible context.
Variant do_call_method(std::string const &method_name, VariantMap const &params) override
std::vector< std::shared_ptr< VTKHandle > > m_vtk_writers
auto find_vtk(std::shared_ptr< VTKHandle > const &vtk) const
virtual void make_instance(VariantMap const &params)=0
virtual ::LatticeModel::units_map get_lattice_to_md_units_conversion() 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_vector_of_variants(std::vector< T > const &v)
Definition Variant.hpp:148
Recursive variant implementation.
Definition Variant.hpp:84