ESPResSo
Extensible Simulation Package for Research on Soft Matter Systems
Loading...
Searching...
No Matches
LatticeModel.cpp
Go to the documentation of this file.
1/*
2 * Copyright (C) 2020-2026 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
22
23#include <blockforest/StructuredBlockForest.h>
24#include <field/vtk/VTKWriter.h>
25#include <vtk/VTKOutput.h>
26
27#include <memory>
28#include <string>
29
30std::shared_ptr<VTKHandle>
33 std::string const &identifier,
34 std::string const &base_folder,
35 std::string const &prefix, bool force_pvtu) {
36
37 using walberla::uint_c;
38
39 // VTKOutput object must be unique
40 auto const vtk_uid = base_folder + "/" + identifier;
41 if (m_vtk_auto.contains(vtk_uid) or m_vtk_manual.contains(vtk_uid)) {
42 throw vtk_runtime_error(vtk_uid, "already exists");
43 }
44
45 // instantiate VTKOutput object
46 auto const &blocks = get_lattice().get_blocks();
47 auto const write_freq = (delta_N) ? static_cast<unsigned int>(delta_N) : 1u;
48 auto vtk_obj = walberla::vtk::createVTKOutput_BlockData(
50 base_folder, prefix, true, true, true, true, uint_c(initial_count));
51
52 // add filters
54
55 // add writers
57
58 auto vtk_handle = std::make_shared<VTKHandle>(vtk_obj, initial_count, true);
59 if (delta_N) {
61 } else {
63 }
64 return vtk_handle;
65}
66
67void LatticeModel::write_vtk(std::string const &vtk_uid) {
68 if (m_vtk_auto.contains(vtk_uid)) {
69 throw vtk_runtime_error(vtk_uid, "is an automatic observable");
70 }
71 if (not m_vtk_manual.contains(vtk_uid)) {
72 throw vtk_runtime_error(vtk_uid, "doesn't exist");
73 }
75 walberla::vtk::writeFiles(vtk_handle->ptr)();
76 vtk_handle->execution_count++;
77}
78
79void LatticeModel::switch_vtk(std::string const &vtk_uid, bool status) {
80 if (m_vtk_manual.contains(vtk_uid)) {
81 throw vtk_runtime_error(vtk_uid, "is a manual observable");
82 }
83 if (not m_vtk_auto.contains(vtk_uid)) {
84 throw vtk_runtime_error(vtk_uid, "doesn't exist");
85 }
86 m_vtk_auto[vtk_uid]->enabled = status;
87}
std::map< std::string, std::shared_ptr< VTKHandle > > m_vtk_auto
VTK writers that are executed automatically.
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.
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, bool force_pvtu)
Create a VTK observable.
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.
cudaStream_t stream[1]
CUDA streams for parallel computing on CPU and GPU.