ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
walberla_bridge/include/walberla_bridge/LatticeModel.hpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2019-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
24
25#include <map>
26#include <memory>
27#include <string>
28#include <unordered_map>
29
30/** @brief Abstract representation of a lattice-based model. */
32public:
33 using units_map = std::unordered_map<std::string, double>;
34
35protected:
36 /** VTK writers that are executed automatically */
37 std::map<std::string, std::shared_ptr<VTKHandle>> m_vtk_auto;
38 /** VTK writers that are executed manually */
39 std::map<std::string, std::shared_ptr<VTKHandle>> m_vtk_manual;
40
41 /** Register VTK writers. Use the multi-piece uniform grid format. */
42 virtual void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj,
43 units_map const &units_conversion,
44 int flag_observables) = 0;
45
46 virtual void
47 register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj) = 0;
48
49 virtual void integrate_vtk_writers() = 0;
50
51public:
52 virtual ~LatticeModel() = default;
53
54 /** @brief Get the underlying lattice. */
55 virtual LatticeWalberla const &get_lattice() const noexcept = 0;
56
57 /** @brief Create a VTK observable.
58 *
59 * @param delta_N Write frequency, if 0 write a single frame,
60 * otherwise add a callback to write every
61 * @p delta_N EK steps to a new file
62 * @param initial_count Initial execution count
63 * @param flag_observables Which observables to measure (OR'ing of
64 * @ref OutputVTK or @ref EKOutputVTK values)
65 * @param units_conversion Lattice-to-MD units conversion
66 * @param identifier Name of the VTK dataset
67 * @param base_folder Path to the VTK folder
68 * @param prefix Prefix of the VTK files
69 */
70 std::shared_ptr<VTKHandle>
71 create_vtk(int delta_N, int initial_count, int flag_observables,
72 units_map const &units_conversion, std::string const &identifier,
73 std::string const &base_folder, std::string const &prefix);
74
75 /** @brief Write a VTK observable to disk.
76 *
77 * @param vtk_uid Name of the VTK object
78 */
79 void write_vtk(std::string const &vtk_uid);
80
81 /** @brief Toggle a VTK observable on/off.
82 *
83 * @param vtk_uid Name of the VTK object
84 * @param status @c true to switch on, @c false to switch off
85 */
86 void switch_vtk(std::string const &vtk_uid, bool status);
87};
Abstract representation of a lattice-based model.
std::map< std::string, std::shared_ptr< VTKHandle > > m_vtk_auto
VTK writers that are executed automatically.
std::shared_ptr< VTKHandle > create_vtk(int delta_N, int initial_count, int flag_observables, units_map const &units_conversion, std::string const &identifier, std::string const &base_folder, std::string const &prefix)
Create a VTK observable.
virtual void register_vtk_field_writers(walberla::vtk::VTKOutput &vtk_obj, units_map const &units_conversion, int flag_observables)=0
Register VTK writers.
virtual void register_vtk_field_filters(walberla::vtk::VTKOutput &vtk_obj)=0
std::map< std::string, std::shared_ptr< VTKHandle > > m_vtk_manual
VTK writers that are executed manually.
virtual void integrate_vtk_writers()=0
virtual ~LatticeModel()=default
std::unordered_map< std::string, double > units_map
void switch_vtk(std::string const &vtk_uid, bool status)
Toggle a VTK observable on/off.
void write_vtk(std::string const &vtk_uid)
Write a VTK observable to disk.
virtual LatticeWalberla const & get_lattice() const noexcept=0
Get the underlying lattice.
Class that runs and controls the BlockForest in waLBerla.