OGS
IntegrationPointWriter.cpp
Go to the documentation of this file.
1
12
13#include <nlohmann/json.hpp>
14
15#include "MeshLib/Mesh.h"
17
18using nlohmann::json;
19
25{
26 auto const& ip_values = writer.values(/*t, x, dof_table*/);
27 assert(ip_values.size() == mesh.getNumberOfElements());
28
29 // create field data and fill it with nodal values, and an offsets cell
30 // array indicating where the cell's integration point data starts.
33 writer.numberOfComponents());
34 field_data.clear();
35
36 for (const auto& element_ip_values : ip_values)
37 {
38 std::copy(element_ip_values.begin(), element_ip_values.end(),
39 std::back_inserter(field_data));
40 }
41
42 return {writer.name(), writer.numberOfComponents(),
43 writer.integrationOrder()};
44}
45
49 MeshLib::Mesh& mesh,
50 std::vector<MeshLib::IntegrationPointMetaData> const& meta_data)
51{
52 json json_meta_data;
53 json_meta_data["integration_point_arrays"] = json::array();
54
55 for (auto const& md : meta_data)
56 {
57 json_meta_data["integration_point_arrays"].push_back(
58 {{"name", md.name},
59 {"number_of_components", md.n_components},
60 {"integration_order", md.integration_order}});
61 }
62
63 // Store the field data.
64 std::string const json_string = json_meta_data.dump();
66 mesh, "IntegrationPointMetaData",
68 dictionary.clear();
69 std::copy(json_string.begin(), json_string.end(),
70 std::back_inserter(dictionary));
71}
72
76 json const& meta_data, std::string const& name)
77{
78 auto const& ip_meta_data = meta_data["integration_point_arrays"];
79 if (auto const it =
80 find_if(cbegin(ip_meta_data), cend(ip_meta_data),
81 [&name](auto const& md) { return md["name"] == name; });
82 it != cend(ip_meta_data))
83 {
84 return {name, (*it)["number_of_components"],
85 (*it)["integration_order"]};
86 }
87 OGS_FATAL("No integration point meta data with name '{:s}' found.", name);
88}
89
90namespace MeshLib
91{
93 MeshLib::Mesh& mesh,
94 std::vector<std::unique_ptr<IntegrationPointWriter>> const&
95 integration_point_writer)
96{
97 std::vector<IntegrationPointMetaData> meta_data;
98 meta_data.reserve(size(integration_point_writer));
99 transform(cbegin(integration_point_writer), cend(integration_point_writer),
100 back_inserter(meta_data),
101 [&](auto const& ip_writer)
102 { return addIntegrationPointData(mesh, *ip_writer); });
103 if (!meta_data.empty())
104 {
105 addIntegrationPointMetaData(mesh, meta_data);
106 }
107}
108
110 MeshLib::Properties const& properties, std::string const& name)
111{
112 if (!properties.existsPropertyVector<char>("IntegrationPointMetaData"))
113 {
114 OGS_FATAL(
115 "Integration point data '{:s}' is present in the vtk field data "
116 "but the required 'IntegrationPointMetaData' array is not "
117 "available.",
118 name);
119 }
120 auto const& mesh_property_ip_meta_data =
121 *properties.template getPropertyVector<char>(
122 "IntegrationPointMetaData");
123
124 if (mesh_property_ip_meta_data.getMeshItemType() !=
126 {
127 OGS_FATAL("IntegrationPointMetaData array must be field data.");
128 }
129
130 // Find the current integration point data entry and extract the
131 // meta data.
132 auto const ip_meta_data = extractIntegrationPointMetaData(
133 json::parse(mesh_property_ip_meta_data.begin(),
134 mesh_property_ip_meta_data.end()),
135 name);
136
137 return ip_meta_data;
138}
139} // namespace MeshLib
#define OGS_FATAL(...)
Definition Error.h:26
static MeshLib::IntegrationPointMetaData addIntegrationPointData(MeshLib::Mesh &mesh, MeshLib::IntegrationPointWriter const &writer)
static MeshLib::IntegrationPointMetaData extractIntegrationPointMetaData(json const &meta_data, std::string const &name)
static void addIntegrationPointMetaData(MeshLib::Mesh &mesh, std::vector< MeshLib::IntegrationPointMetaData > const &meta_data)
Definition of the Mesh class.
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:97
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
Definition Properties.h:36
bool existsPropertyVector(std::string_view name) const
Definition Properties.h:74
PropertyVector< T > * getOrCreateMeshProperty(Mesh &mesh, std::string const &property_name, MeshItemType const item_type, int const number_of_components)
IntegrationPointMetaData getIntegrationPointMetaData(MeshLib::Properties const &properties, std::string const &name)
void addIntegrationPointDataToMesh(MeshLib::Mesh &mesh, std::vector< std::unique_ptr< IntegrationPointWriter > > const &integration_point_writer)
std::vector< std::vector< double > > values() const