OGS
OGSMesh.cpp
Go to the documentation of this file.
1
13#include "OGSMesh.h"
14
15#include <algorithm>
16#include <pybind11/eigen.h>
17#include <pybind11/pybind11.h>
18#include <pybind11/stl.h>
19#include <spdlog/spdlog.h>
20
21#include <numeric>
22#include <range/v3/algorithm/copy.hpp>
23#include <range/v3/numeric.hpp>
24#include <range/v3/range/conversion.hpp>
25#include <range/v3/view/enumerate.hpp>
26#include <range/v3/view/indirect.hpp>
27#include <range/v3/view/join.hpp>
28#include <range/v3/view/map.hpp>
29#include <range/v3/view/transform.hpp>
30#include <span>
31#include <vector>
32
34#include "BaseLib/Logging.h"
35#include "InfoLib/GitInfo.h"
37#include "MeshLib/Mesh.h"
38#include "MeshLib/Node.h"
40#include "MeshLib/VtkOGSEnum.h"
41
42OGSMesh::OGSMesh(MeshLib::Mesh& mesh) : _mesh(mesh) {}
43
44std::vector<double> OGSMesh::getPointCoordinates() const
45{
46 return ranges::to<std::vector>(_mesh.getNodes() | MeshLib::views::coords |
47 ranges::views::join);
48}
49
50std::pair<std::vector<int>, std::vector<int>> OGSMesh::getCells() const
51{
52 auto const& elements = _mesh.getElements();
53 std::vector<int> cells;
54 std::vector<int> cell_types;
55 for (auto const* element : elements)
56 {
57 cells.push_back(static_cast<int>(element->getNumberOfNodes()));
58 ranges::copy(element->nodes() | MeshLib::views::ids,
59 std::back_inserter(cells));
60 cell_types.push_back(OGSToVtkCellType(element->getCellType()));
61 }
62 return {cells, cell_types};
63}
64
65void OGSMesh::setPointDataArray(std::string const& name,
66 std::vector<double> const& values,
67 std::size_t const number_of_components)
68{
69 auto* pv = MeshLib::getOrCreateMeshProperty<double>(
70 _mesh, name, MeshLib::MeshItemType::Node, number_of_components);
71 if (pv == nullptr)
72 {
73 OGS_FATAL("Couldn't access cell/element property '{}'.", name);
74 }
75 if (pv->size() != values.size())
76 {
78 "OGSMesh::setPointDataArray: size mismatch: property vector has "
79 "size '{}', while the number of values is '{}'.",
80 pv->size(), values.size());
81 }
82 std::copy(values.begin(), values.end(), pv->data());
83}
84
85std::vector<double> OGSMesh::getPointDataArray(
86 std::string const& name, std::size_t const number_of_components) const
87{
88 auto const* pv = MeshLib::getOrCreateMeshProperty<double>(
89 _mesh, name, MeshLib::MeshItemType::Node, number_of_components);
90 if (pv == nullptr)
91 {
92 OGS_FATAL("Couldn't access point/node property '{}'.", name);
93 }
94 return ranges::to<std::vector>(std::span(pv->data(), pv->size()));
95}
96
97void OGSMesh::setCellDataArray(std::string const& name,
98 std::vector<double> const& values,
99 std::size_t const number_of_components)
100{
101 auto* pv = MeshLib::getOrCreateMeshProperty<double>(
102 _mesh, name, MeshLib::MeshItemType::Cell, number_of_components);
103 if (pv == nullptr)
104 {
105 OGS_FATAL("Couldn't access cell/element property '{}'.", name);
106 }
107 if (pv->size() != values.size())
108 {
109 OGS_FATAL(
110 "OGSMesh::setCellDataArray: size mismatch: property vector has "
111 "size '{}', while the number of values is '{}'.",
112 pv->size(), values.size());
113 }
114 std::copy(values.begin(), values.end(), pv->data());
115}
116
117std::vector<double> OGSMesh::getCellDataArray(
118 std::string const& name, std::size_t const number_of_components) const
119{
120 auto const* pv = MeshLib::getOrCreateMeshProperty<double>(
121 _mesh, name, MeshLib::MeshItemType::Cell, number_of_components);
122 if (pv == nullptr)
123 {
124 OGS_FATAL("Couldn't access cell/element property '{}'.", name);
125 }
126 return ranges::to<std::vector>(std::span(pv->data(), pv->size()));
127}
Definition of the Element class.
#define OGS_FATAL(...)
Definition Error.h:26
Git information.
Definition of the Mesh class.
Definition of the Node class.
int OGSToVtkCellType(MeshLib::CellType ogs)
std::vector< Node * > const & getNodes() const
Get the nodes-vector for the mesh.
Definition Mesh.h:106
std::vector< Element * > const & getElements() const
Get the element-vector for the mesh.
Definition Mesh.h:109
void setPointDataArray(std::string const &name, std::vector< double > const &values, std::size_t const number_of_components)
Definition OGSMesh.cpp:65
std::vector< double > getPointDataArray(std::string const &name, std::size_t const number_of_components=1) const
Definition OGSMesh.cpp:85
void setCellDataArray(std::string const &name, std::vector< double > const &values, std::size_t const number_of_components)
Definition OGSMesh.cpp:97
MeshLib::Mesh & _mesh
Definition OGSMesh.h:46
std::vector< double > getPointCoordinates() const
Definition OGSMesh.cpp:44
OGSMesh(MeshLib::Mesh &mesh)
Definition OGSMesh.cpp:42
std::pair< std::vector< int >, std::vector< int > > getCells() const
Definition OGSMesh.cpp:50
std::vector< double > getCellDataArray(std::string const &name, std::size_t const number_of_components) const
Definition OGSMesh.cpp:117
constexpr ranges::views::view_closure ids
For an element of a range view return its id.
Definition Mesh.h:225
constexpr ranges::views::view_closure coords
Definition Mesh.h:232