OGS
MeshToolsLib::RasterToMesh Class Reference

Detailed Description

Converts raster data into an OGS mesh.

Definition at line 28 of file RasterToMesh.h.

#include <RasterToMesh.h>

Static Public Member Functions

static std::unique_ptr< MeshLib::Meshconvert (GeoLib::Raster const &raster, MeshLib::MeshElemType elem_type, MeshLib::UseIntensityAs intensity_type, std::string const &array_name="Colour")
static std::unique_ptr< MeshLib::Meshconvert (vtkImageData *img, const double origin[3], const double scalingFactor, MeshLib::MeshElemType elem_type, MeshLib::UseIntensityAs intensity_type, std::string const &array_name="Colour")
static std::unique_ptr< MeshLib::Meshconvert (const double *const img, GeoLib::RasterHeader const &header, MeshLib::MeshElemType elem_type, MeshLib::UseIntensityAs intensity_type, std::string const &array_name="Colour")

Static Private Member Functions

template<typename T>
static void fillPropertyVector (MeshLib::PropertyVector< T > &prop_vec, double const *const img, GeoLib::RasterHeader const &header, MeshLib::MeshElemType elem_type)

Member Function Documentation

◆ convert() [1/3]

std::unique_ptr< MeshLib::Mesh > MeshToolsLib::RasterToMesh::convert ( const double *const img,
GeoLib::RasterHeader const & header,
MeshLib::MeshElemType elem_type,
MeshLib::UseIntensityAs intensity_type,
std::string const & array_name = "Colour" )
static

Converts double array with raster values into a mesh.

Parameters
imginput image.
headerraster header information.
elem_typedefines if elements of the new mesh should be triangles or quads (or hexes for 3D).
intensity_typedefines how image intensities are interpreted.
array_namemesh property name, defaults to "Colour" if not given.

Definition at line 37 of file RasterToMesh.cpp.

43{
44 if ((elem_type != MeshLib::MeshElemType::TRIANGLE) &&
45 (elem_type != MeshLib::MeshElemType::QUAD) &&
46 (elem_type != MeshLib::MeshElemType::HEXAHEDRON) &&
47 (elem_type != MeshLib::MeshElemType::PRISM))
48 {
49 ERR("Invalid Mesh Element Type.");
50 return nullptr;
51 }
52
53 if (((elem_type == MeshLib::MeshElemType::TRIANGLE) ||
54 (elem_type == MeshLib::MeshElemType::QUAD)) &&
55 header.n_depth != 1)
56 {
57 ERR("Triangle or Quad elements cannot be used to construct meshes from "
58 "3D rasters.");
59 return nullptr;
60 }
61
62 if (intensity_type == MeshLib::UseIntensityAs::ELEVATION &&
63 ((elem_type == MeshLib::MeshElemType::PRISM) ||
65 {
66 ERR("Elevation mapping can only be performed for 2D meshes.");
67 return nullptr;
68 }
69
70 std::unique_ptr<MeshLib::Mesh> mesh;
71 if (elem_type == MeshLib::MeshElemType::TRIANGLE)
72 {
74 header.n_cols,
75 header.n_rows,
76 header.cell_size,
77 header.origin,
78 "RasterDataMesh"));
79 }
80 else if (elem_type == MeshLib::MeshElemType::QUAD)
81 {
83 header.n_cols,
84 header.n_rows,
85 header.cell_size,
86 header.origin,
87 "RasterDataMesh"));
88 }
89 else if (elem_type == MeshLib::MeshElemType::PRISM)
90 {
92 header.n_cols,
93 header.n_rows,
94 header.n_depth,
95 header.cell_size,
96 header.origin,
97 "RasterDataMesh"));
98 }
99 else if (elem_type == MeshLib::MeshElemType::HEXAHEDRON)
100 {
102 header.n_cols,
103 header.n_rows,
104 header.n_depth,
105 header.cell_size,
106 header.origin,
107 "RasterDataMesh"));
108 }
109
110 std::unique_ptr<MeshLib::Mesh> new_mesh;
111 std::vector<std::size_t> elements_to_remove;
112 if (intensity_type == MeshLib::UseIntensityAs::ELEVATION)
113 {
114 std::vector<MeshLib::Node*> const& nodes(mesh->getNodes());
115 std::vector<MeshLib::Element*> const& elems(mesh->getElements());
116 std::size_t const n_nodes(elems[0]->getNumberOfNodes());
117 bool const double_idx =
118 (elem_type == MeshLib::MeshElemType::TRIANGLE) ||
119 (elem_type == MeshLib::MeshElemType::PRISM);
120 std::size_t const m = (double_idx) ? 2 : 1;
121 for (std::size_t k = 0; k < header.n_depth; k++)
122 {
123 std::size_t const layer_idx = (k * header.n_rows * header.n_cols);
124 for (std::size_t i = 0; i < header.n_cols; i++)
125 {
126 std::size_t const idx(i * header.n_rows + layer_idx);
127 for (std::size_t j = 0; j < header.n_rows; j++)
128 {
129 double const val(img[idx + j]);
130 if (val == header.no_data)
131 {
132 elements_to_remove.push_back(m * (idx + j));
133 if (double_idx)
134 {
135 elements_to_remove.push_back(m * (idx + j) + 1);
136 }
137 continue;
138 }
139 for (std::size_t n = 0; n < n_nodes; ++n)
140 {
141 (*(nodes[getNodeIndex(*elems[m * (idx + j)], n)]))[2] =
142 val;
143 if (double_idx)
144 {
145 (*(nodes[getNodeIndex(*elems[m * (idx + j) + 1],
146 n)]))[2] = val;
147 }
148 }
149 }
150 }
151 }
152 }
153 else
154 {
155 MeshLib::Properties& properties = mesh->getProperties();
156 MeshLib::ElementSearch ex(*mesh);
157 if (array_name == "MaterialIDs")
158 {
159 auto* const prop_vec = properties.createNewPropertyVector<int>(
160 array_name, MeshLib::MeshItemType::Cell, 1);
161 fillPropertyVector<int>(*prop_vec, img, header, elem_type);
162 ex.searchByPropertyValue<int>(array_name,
163 static_cast<int>(header.no_data));
164 }
165 else
166 {
167 auto* const prop_vec = properties.createNewPropertyVector<double>(
168 array_name, MeshLib::MeshItemType::Cell, 1);
169 fillPropertyVector<double>(*prop_vec, img, header, elem_type);
170 ex.searchByPropertyValue<double>(array_name, header.no_data);
171 }
172 elements_to_remove = ex.getSearchedElementIDs();
173 }
174 if (!elements_to_remove.empty())
175 {
176 new_mesh.reset(MeshToolsLib::removeElements(
177 *mesh, elements_to_remove, mesh->getName()));
178 }
179 else
180 {
181 new_mesh = std::move(mesh);
182 }
183
184 if (intensity_type == MeshLib::UseIntensityAs::NONE)
185 {
186 new_mesh->getProperties().removePropertyVector(array_name);
187 }
188
189 return new_mesh;
190}
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
PropertyVector< T > * createNewPropertyVector(std::string_view name, MeshItemType mesh_item_type, std::size_t n_components=1)
static void fillPropertyVector(MeshLib::PropertyVector< T > &prop_vec, double const *const img, GeoLib::RasterHeader const &header, MeshLib::MeshElemType elem_type)
std::size_t getNodeIndex(Element const &element, unsigned const idx)
Definition Element.cpp:226
MeshLib::Mesh * generateRegularQuadMesh(const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::Mesh * generateRegularPrismMesh(const double x_length, const double y_length, const double z_length, const std::size_t x_subdivision, const std::size_t y_subdivision, const std::size_t z_subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::Mesh * generateRegularTriMesh(const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::Mesh * generateRegularHexMesh(const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, const BaseLib::ISubdivision &div_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::Mesh * removeElements(const MeshLib::Mesh &mesh, const std::vector< std::size_t > &removed_element_ids, const std::string &new_mesh_name)

References MeshLib::Cell, GeoLib::RasterHeader::cell_size, MeshLib::Properties::createNewPropertyVector(), MeshLib::ELEVATION, ERR(), fillPropertyVector(), MeshToolsLib::MeshGenerator::generateRegularHexMesh(), MeshToolsLib::MeshGenerator::generateRegularPrismMesh(), MeshToolsLib::MeshGenerator::generateRegularQuadMesh(), MeshToolsLib::MeshGenerator::generateRegularTriMesh(), MeshLib::ElementSearch::getSearchedElementIDs(), MeshLib::HEXAHEDRON, GeoLib::RasterHeader::n_cols, GeoLib::RasterHeader::n_depth, GeoLib::RasterHeader::n_rows, GeoLib::RasterHeader::no_data, MeshLib::NONE, GeoLib::RasterHeader::origin, MeshLib::PRISM, MeshLib::QUAD, MeshToolsLib::removeElements(), MeshLib::ElementSearch::searchByPropertyValue(), and MeshLib::TRIANGLE.

◆ convert() [2/3]

std::unique_ptr< MeshLib::Mesh > MeshToolsLib::RasterToMesh::convert ( GeoLib::Raster const & raster,
MeshLib::MeshElemType elem_type,
MeshLib::UseIntensityAs intensity_type,
std::string const & array_name = "Colour" )
static

Converts greyscale raster into a mesh.

Parameters
rasterinput raster.
elem_typedefines if elements of the new mesh should be triangles or quads (or hexes for 3D).
intensity_typedefines how image intensities are interpreted.
array_namemesh property name, defaults to "Colour" if not given.

Definition at line 24 of file RasterToMesh.cpp.

29{
30 return convert(raster.data(),
31 raster.getHeader(),
32 elem_type,
33 intensity_type,
34 array_name);
35}
static std::unique_ptr< MeshLib::Mesh > convert(GeoLib::Raster const &raster, MeshLib::MeshElemType elem_type, MeshLib::UseIntensityAs intensity_type, std::string const &array_name="Colour")

References convert(), GeoLib::Raster::data(), and GeoLib::Raster::getHeader().

Referenced by convert(), convert(), NetCdfConfigureDialog::createDataObject(), main(), VtkVisPipelineView::showImageToMeshConversionDialog(), and writeDataToMesh().

◆ convert() [3/3]

std::unique_ptr< MeshLib::Mesh > MeshToolsLib::RasterToMesh::convert ( vtkImageData * img,
const double origin[3],
const double scalingFactor,
MeshLib::MeshElemType elem_type,
MeshLib::UseIntensityAs intensity_type,
std::string const & array_name = "Colour" )
static

Converts a vtkImageData into a mesh.

Parameters
imginput image.
origincoordinates of image's origin, lower left corner.
scalingFactoredge length of each pixel
elem_typedefines if elements of the new mesh should be triangles or quads (or hexes for 3D).
intensity_typedefines how image intensities are interpreted.
array_namemesh property name, defaults to "Colour" if not given.

◆ fillPropertyVector()

template<typename T>
void MeshToolsLib::RasterToMesh::fillPropertyVector ( MeshLib::PropertyVector< T > & prop_vec,
double const *const img,
GeoLib::RasterHeader const & header,
MeshLib::MeshElemType elem_type )
inlinestaticprivate

Definition at line 84 of file RasterToMesh.h.

88 {
89 for (std::size_t k = 0; k < header.n_depth; k++)
90 {
91 std::size_t const layer_idx = (k * header.n_rows * header.n_cols);
92 for (std::size_t i = 0; i < header.n_cols; i++)
93 {
94 std::size_t idx(i * header.n_rows + layer_idx);
95 for (std::size_t j = 0; j < header.n_rows; j++)
96 {
97 auto val(static_cast<T>(img[idx + j]));
98 prop_vec.push_back(val);
99 if (elem_type == MeshLib::MeshElemType::TRIANGLE ||
100 elem_type == MeshLib::MeshElemType::PRISM)
101 {
102 prop_vec.push_back(val); // because each pixel is
103 // represented by two cells
104 }
105 }
106 }
107 }
108 }
constexpr void push_back(const PROP_VAL_TYPE &value)

References GeoLib::RasterHeader::n_cols, GeoLib::RasterHeader::n_depth, GeoLib::RasterHeader::n_rows, MeshLib::PRISM, MeshLib::PropertyVector< PROP_VAL_TYPE >::push_back(), and MeshLib::TRIANGLE.

Referenced by convert().


The documentation for this class was generated from the following files: