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
28
29#include <algorithm>
30#include <memory>
31#include <stdexcept>
32#include <string>
33#include <vector>
34
36
37template <class Method, class VTKHandle>
38class LatticeModel : public AutoParameters<LatticeModel<Method, VTKHandle>> {
39protected:
40 std::shared_ptr<LatticeWalberla> m_lattice;
41 std::shared_ptr<Method> m_instance;
42 std::vector<std::shared_ptr<VTKHandle>> m_vtk_writers;
43
44 virtual ::LatticeModel::units_map
46 virtual void make_instance(VariantMap const &params) = 0;
47
48 auto find_vtk(std::shared_ptr<VTKHandle> const &vtk) const {
49 return std::find(m_vtk_writers.begin(), m_vtk_writers.end(), vtk);
50 }
51
55
56public:
57 Variant do_call_method(std::string const &method_name,
58 VariantMap const &params) override {
59 if (method_name == "add_vtk_writer") {
60 auto vtk = get_value<std::shared_ptr<VTKHandle>>(params, "vtk");
61 auto const needle = find_vtk(vtk);
63 if (needle != m_vtk_writers.end()) {
64 throw std::runtime_error(
65 "VTK object is already attached to this lattice");
66 }
67 vtk->attach_to_lattice(m_instance, get_latice_to_md_units_conversion());
68 m_vtk_writers.emplace_back(vtk);
69 });
70 return {};
71 }
72 if (method_name == "remove_vtk_writer") {
73 auto const vtk = get_value<std::shared_ptr<VTKHandle>>(params, "vtk");
74 auto const needle = find_vtk(vtk);
76 if (needle == m_vtk_writers.end()) {
77 throw std::runtime_error(
78 "VTK object is not attached to this lattice");
79 }
80 vtk->detach_from_lattice();
81 });
82 m_vtk_writers.erase(needle);
83 return {};
84 }
85 if (method_name == "clear_vtk_writers") {
86 for (auto const &vtk : m_vtk_writers) {
87 vtk->detach_from_lattice();
88 }
89 m_vtk_writers.clear();
90 }
91 return {};
92 }
93};
94
95} // 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.
virtual ::LatticeModel::units_map get_latice_to_md_units_conversion() const =0
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
std::unordered_map< std::string, Variant > VariantMap
Definition Variant.hpp:82
auto make_vector_of_variants(std::vector< T > const &v)
Definition Variant.hpp:101
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:80
static SteepestDescentParameters params
Currently active steepest descent instance.