OGS
IntegrationPointWriter.cpp
Go to the documentation of this file.
1 
11 #include "IntegrationPointWriter.h"
12 
13 #include <nlohmann/json.hpp>
14 
15 #include "MeshLib/Mesh.h"
16 
17 using nlohmann::json;
18 
24 {
25  auto const& ip_values = writer.values(/*t, x, dof_table*/);
26  assert(ip_values.size() == mesh.getNumberOfElements());
27 
28  // create field data and fill it with nodal values, and an offsets cell
29  // array indicating where the cell's integration point data starts.
30  auto& field_data = *MeshLib::getOrCreateMeshProperty<double>(
32  writer.numberOfComponents());
33  field_data.clear();
34 
35  for (const auto& element_ip_values : ip_values)
36  {
37  std::copy(element_ip_values.begin(), element_ip_values.end(),
38  std::back_inserter(field_data));
39  }
40 
41  return {writer.name(), writer.numberOfComponents(),
42  writer.integrationOrder()};
43 }
44 
48  MeshLib::Mesh& mesh,
49  std::vector<ProcessLib::IntegrationPointMetaData> const& meta_data)
50 {
51  json json_meta_data;
52  json_meta_data["integration_point_arrays"] = json::array();
53 
54  for (auto const& md : meta_data)
55  {
56  json_meta_data["integration_point_arrays"].push_back(
57  {{"name", md.name},
58  {"number_of_components", md.n_components},
59  {"integration_order", md.integration_order}});
60  }
61 
62  // Store the field data.
63  std::string const json_string = json_meta_data.dump();
64  auto& dictionary = *MeshLib::getOrCreateMeshProperty<char>(
65  mesh, "IntegrationPointMetaData",
67  dictionary.clear();
68  std::copy(json_string.begin(), json_string.end(),
69  std::back_inserter(dictionary));
70 }
71 
75  json const& meta_data, std::string const& name)
76 {
77  auto const& ip_meta_data = meta_data["integration_point_arrays"];
78  if (auto const it =
79  find_if(cbegin(ip_meta_data), cend(ip_meta_data),
80  [&name](auto const& md) { return md["name"] == name; });
81  it != cend(ip_meta_data))
82  {
83  return {name, (*it)["number_of_components"],
84  (*it)["integration_order"]};
85  }
86  OGS_FATAL("No integration point meta data with name '{:s}' found.", name);
87 }
88 
89 namespace ProcessLib
90 {
92  MeshLib::Mesh& mesh,
93  std::vector<std::unique_ptr<IntegrationPointWriter>> const&
94  integration_point_writer)
95 {
96  std::vector<IntegrationPointMetaData> meta_data;
97  meta_data.reserve(size(integration_point_writer));
98  transform(cbegin(integration_point_writer), cend(integration_point_writer),
99  back_inserter(meta_data),
100  [&](auto const& ip_writer)
101  { return addIntegrationPointData(mesh, *ip_writer); });
102  if (!meta_data.empty())
103  {
104  addIntegrationPointMetaData(mesh, meta_data);
105  }
106 }
107 
109  std::string const& name)
110 {
111  if (!mesh.getProperties().existsPropertyVector<char>(
112  "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  *mesh.getProperties().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 ProcessLib
#define OGS_FATAL(...)
Definition: Error.h:26
static ProcessLib::IntegrationPointMetaData extractIntegrationPointMetaData(json const &meta_data, std::string const &name)
static ProcessLib::IntegrationPointMetaData addIntegrationPointData(MeshLib::Mesh &mesh, ProcessLib::IntegrationPointWriter const &writer)
static void addIntegrationPointMetaData(MeshLib::Mesh &mesh, std::vector< ProcessLib::IntegrationPointMetaData > const &meta_data)
Definition of the Mesh class.
Properties & getProperties()
Definition: Mesh.h:123
std::size_t getNumberOfElements() const
Get the number of elements.
Definition: Mesh.h:86
bool existsPropertyVector(std::string const &name) const
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37
IntegrationPointMetaData getIntegrationPointMetaData(MeshLib::Mesh const &mesh, std::string const &name)
void addIntegrationPointWriter(MeshLib::Mesh &mesh, std::vector< std::unique_ptr< IntegrationPointWriter >> const &integration_point_writer)
std::vector< std::vector< double > > values() const