OGS
IntegrationPointWriter.cpp File Reference

Detailed Description

Definition in file IntegrationPointWriter.cpp.

#include "IntegrationPointWriter.h"
#include <nlohmann/json.hpp>
#include "MeshLib/Mesh.h"
#include "MeshLib/Utils/getOrCreateMeshProperty.h"
Include dependency graph for IntegrationPointWriter.cpp:

Go to the source code of this file.

Namespaces

namespace  MeshLib
 

Functions

static MeshLib::IntegrationPointMetaData addIntegrationPointData (MeshLib::Mesh &mesh, MeshLib::IntegrationPointWriter const &writer)
 
static void addIntegrationPointMetaData (MeshLib::Mesh &mesh, std::vector< MeshLib::IntegrationPointMetaData > const &meta_data)
 
static MeshLib::IntegrationPointMetaData extractIntegrationPointMetaData (json const &meta_data, std::string const &name)
 
void MeshLib::addIntegrationPointDataToMesh (MeshLib::Mesh &mesh, std::vector< std::unique_ptr< IntegrationPointWriter > > const &integration_point_writer)
 
IntegrationPointMetaData MeshLib::getIntegrationPointMetaData (MeshLib::Properties const &properties, std::string const &name)
 

Function Documentation

◆ addIntegrationPointData()

static MeshLib::IntegrationPointMetaData addIntegrationPointData ( MeshLib::Mesh & mesh,
MeshLib::IntegrationPointWriter const & writer )
static

Adds the integration point data and creates meta data for it.

Returns meta data for the written integration point data.

Definition at line 23 of file IntegrationPointWriter.cpp.

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.
31 auto& field_data = *MeshLib::getOrCreateMeshProperty<double>(
32 mesh, writer.name(), MeshLib::MeshItemType::IntegrationPoint,
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}
std::size_t getNumberOfElements() const
Get the number of elements.
Definition Mesh.h:97

References MeshLib::Mesh::getNumberOfElements(), MeshLib::IntegrationPointWriter::integrationOrder(), MeshLib::IntegrationPoint, MeshLib::IntegrationPointWriter::name(), MeshLib::IntegrationPointWriter::numberOfComponents(), and MeshLib::IntegrationPointWriter::values().

Referenced by MeshLib::addIntegrationPointDataToMesh().

◆ addIntegrationPointMetaData()

static void addIntegrationPointMetaData ( MeshLib::Mesh & mesh,
std::vector< MeshLib::IntegrationPointMetaData > const & meta_data )
static

Adds integration point meta data as char mesh property encoded in JSON format, which is then stored as VTK's field data.

Definition at line 48 of file IntegrationPointWriter.cpp.

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();
65 auto& dictionary = *MeshLib::getOrCreateMeshProperty<char>(
66 mesh, "IntegrationPointMetaData",
68 dictionary.clear();
69 std::copy(json_string.begin(), json_string.end(),
70 std::back_inserter(dictionary));
71}

References MeshLib::IntegrationPoint.

Referenced by MeshLib::addIntegrationPointDataToMesh().

◆ extractIntegrationPointMetaData()

static MeshLib::IntegrationPointMetaData extractIntegrationPointMetaData ( json const & meta_data,
std::string const & name )
static

For the given json object and the name extract integration point meta data, or fail if no meta data was found for the given name.

Definition at line 75 of file IntegrationPointWriter.cpp.

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}
#define OGS_FATAL(...)
Definition Error.h:26

References OGS_FATAL.

Referenced by MeshLib::getIntegrationPointMetaData().