CortidQCT  1.2.2.52
Mesh.h
Go to the documentation of this file.
1 
12 #pragma once
13 
14 #include "BarycentricPoint.h"
15 #include "ColorToLabelMap.h"
16 #include "LabelToColorMap.h"
17 #include "Ray.h"
18 #include "RayMeshIntersection.h"
19 
20 #include <array>
21 #include <stdexcept>
22 #include <string>
23 #include <type_traits>
24 #include <vector>
25 
26 namespace CortidQCT {
27 
35 template <class T> class Mesh {
36 
37  static_assert(std::is_floating_point<T>::value,
38  "Scalar must be a floating point type");
39 
40 public:
42  using Scalar = T;
43 
45  using Index = std::ptrdiff_t;
46 
48  using Size = std::size_t;
49 
51  using Label = unsigned int;
52 
53 private:
55  using VertexData = std::vector<Scalar>;
57  using IndexData = std::vector<Index>;
59  using LabelData = std::vector<Label>;
61  using NormalData = VertexData;
62 
63 public:
66 
68  inline Mesh() noexcept(noexcept(VertexData()) && noexcept(IndexData()) &&
69  noexcept(LabelData()) && noexcept(VertexData())) {}
70 
72  inline Mesh(std::size_t nVertices, std::size_t nTriangles) noexcept(
73  noexcept(VertexData(3 * nVertices, T{0})) &&
74  noexcept(IndexData(3 * nTriangles, Index{0})) &&
75  noexcept(LabelData(nVertices, Label{0})) &&
76  noexcept(NormalData(3 * nVertices, T{0})))
77  : vertexData_(3 * nVertices, T{0}), indexData_(3 * nTriangles, Index{0}),
78  labelData_(nVertices, Label{0}), normalData_(3 * nVertices, T{0}) {}
79 
80  // @}
81 
84 
86  inline Size vertexCount() const noexcept { return vertexData_.size() / 3; }
87 
89  inline Size triangleCount() const noexcept { return indexData_.size() / 3; }
90 
94  inline bool isEmpty() const noexcept {
95  return triangleCount() == 0 || vertexCount() == 0;
96  }
98 
101 
118  Mesh &loadFromFile(std::string const &meshFilename,
119  std::string const &labelFilename);
120 
138  Mesh &loadFromFile(std::string const &meshFilename,
139  ColorToLabelMap<Label, double> const &colorMap =
140  ColorToLabelMaps::defaultMap<Label, double>);
141 
157  void writeToFile(std::string const &meshFilename,
158  std::string const &labelFilename) const;
159 
179  void writeToFile(std::string const &meshFilename,
180  LabelToColorMap<double, Label> const &labelMap =
181  LabelToColorMaps::defaultMap<double, Label>) const;
182 
184 
197  std::array<T, 3>
199 
213  template <class InputIterator, class OutputIterator>
214  void cartesianRepresentation(InputIterator begin, InputIterator end,
215  OutputIterator out) const;
216 
224  inline std::vector<std::array<T, 3>> cartesianRepresentation(
225  std::vector<BarycentricPoint<T, Index>> const &points) const {
226  using Cartesian = std::array<T, 3>;
227  std::vector<Cartesian> cartPoints(points.size());
228 
229  auto const out = reinterpret_cast<T *>(cartPoints.data());
230  cartesianRepresentation(points.cbegin(), points.cend(), out);
231 
232  return cartPoints;
233  }
234 
235 #pragma clang diagnostic push
236 #pragma clang diagnostic ignored "-Wdocumentation"
237 
267  template <class PtIter, class AttrIter, class OutputIterator>
268  void barycentricInterpolation(PtIter pointsBegin, PtIter pointsEnd,
269  AttrIter attributesBegin, OutputIterator out,
270  std::size_t attributeDimension = 1) const;
271 #pragma clang diagnostic pop
272 
273 #pragma clang diagnostic push
274 #pragma clang diagnostic ignored "-Wdocumentation"
275 
289  template <class InputIterator, class OutputIterator>
290  void rayIntersections(InputIterator raysBegin, InputIterator raysEnd,
291  OutputIterator intersectionsOut) const;
292 
301 
305  void updatePerVertexNormals();
306 
307 #pragma clang diagnostic pop
308 
310 
322  Mesh<T> &upsample(std::size_t nTimes = 1);
323 
325 
348  template <class F>
349  inline auto withUnsafeVertexPointer(F &&f) const
350  noexcept(noexcept(f(std::declval<const VertexData>().data()))) {
351  return f(vertexData_.data());
352  }
353 
354  template <class F>
355  inline auto withUnsafeVertexPointer(F &&f) noexcept(
356  noexcept(f(std::declval<VertexData>().data()))) {
357  return f(vertexData_.data());
358  }
359 
371  template <class F>
372  inline auto withUnsafeIndexPointer(F &&f) const
373  noexcept(noexcept(f(std::declval<const IndexData>().data()))) {
374  return f(indexData_.data());
375  }
376 
377  template <class F>
378  inline auto withUnsafeIndexPointer(F &&f) noexcept(
379  noexcept(f(std::declval<IndexData>().data()))) {
380  return f(indexData_.data());
381  }
382 
394  template <class F>
395  inline auto withUnsafeLabelPointer(F &&f) const
396  noexcept(noexcept(f(LabelData().data()))) {
397  return f(labelData_.data());
398  }
399 
400  template <class F>
401  inline auto
402  withUnsafeLabelPointer(F &&f) noexcept(noexcept(f(LabelData().data()))) {
403  return f(labelData_.data());
404  }
405 
419  template <class F>
420  inline auto withUnsafeVertexNormalPointer(F &&f) const
421  noexcept(noexcept(f(std::declval<const NormalData>().data()))) {
422  return f(normalData_.data());
423  }
424 
425  template <class F>
426  inline auto withUnsafeVertexNormalPointer(F &&f) noexcept(
427  noexcept(f(std::declval<NormalData>().data()))) {
428  return f(normalData_.data());
429  }
430 
432 private:
434  void ensurePostconditions() const;
435 
437  VertexData vertexData_;
439  IndexData indexData_;
441  LabelData labelData_;
443  NormalData normalData_;
444 };
445 
446 /*************************************
447  * Explicit template instanciations
448  */
449 
450 extern template class Mesh<float>;
451 extern template class Mesh<double>;
452 
453 extern template void Mesh<float>::barycentricInterpolation(
455  BarycentricPoint<float, std::ptrdiff_t> const *, float const *, float *,
456  std::size_t) const;
457 extern template void Mesh<double>::barycentricInterpolation(
459  BarycentricPoint<double, std::ptrdiff_t> const *, double const *, double *,
460  std::size_t) const;
461 
462 extern template void Mesh<float>::cartesianRepresentation(
464  BarycentricPoint<float, std::ptrdiff_t> const *, float *) const;
465 extern template void Mesh<double>::cartesianRepresentation(
467  BarycentricPoint<float, std::ptrdiff_t> const *, double *) const;
468 extern template void Mesh<float>::cartesianRepresentation(
469  std::vector<BarycentricPoint<float, std::ptrdiff_t>>::const_iterator,
470  std::vector<BarycentricPoint<float, std::ptrdiff_t>>::const_iterator,
471  float *) const;
472 extern template void Mesh<double>::cartesianRepresentation(
473  std::vector<BarycentricPoint<double, std::ptrdiff_t>>::const_iterator,
474  std::vector<BarycentricPoint<double, std::ptrdiff_t>>::const_iterator,
475  double *) const;
476 
477 extern template void
480 extern template void Mesh<float>::rayIntersections(
481  Ray<float> const *, Ray<float> const *,
482  std::back_insert_iterator<std::vector<RayMeshIntersection<float>>>) const;
483 extern template void
486 extern template void Mesh<double>::rayIntersections(
487  Ray<double> const *, Ray<double> const *,
488  std::back_insert_iterator<std::vector<RayMeshIntersection<double>>>) const;
489 
490 } // namespace CortidQCT
Name namespace for CortidQCT library.
Definition: CortidQCT.h:23
auto withUnsafeIndexPointer(F &&f) const noexcept(noexcept(f(std::declval< const IndexData >().data())))
Calls the given functional with an unsafe pointer to the raw index storage.
Definition: Mesh.h:372
This finle contains definitions for color to index map types.
std::ptrdiff_t Index
Index type.
Definition: Mesh.h:45
auto withUnsafeVertexPointer(F &&f) const noexcept(noexcept(f(std::declval< const VertexData >().data())))
Calls the given functional with an unsafe pointer to the raw vertex storage.
Definition: Mesh.h:349
void barycentricInterpolation(PtIter pointsBegin, PtIter pointsEnd, AttrIter attributesBegin, OutputIterator out, std::size_t attributeDimension=1) const
Interpolates per-vertex values for points inside a triangle.
Definition: Mesh.cpp:502
void writeToFile(std::string const &meshFilename, std::string const &labelFilename) const
Writes mesh to ASCII file using format auto detection.
Definition: Mesh.cpp:299
A triangle mesh class.
Definition: BarycentricPoint.h:20
std::function< std::array< Scalar, 3 >(Label)> LabelToColorMap
Label to color map.
Definition: LabelToColorMap.h:30
std::array< T, 3 > cartesianRepresentation(BarycentricPoint< T, Index > const &point) const
Returns the cartesian coordinate of a point in barycentric coordinates.
void rayIntersections(InputIterator raysBegin, InputIterator raysEnd, OutputIterator intersectionsOut) const
Computes the intersection of a set of rays with the mesh.
Definition: Mesh.cpp:561
void updatePerVertexNormals()
Re-computes per-vertex normals.
Definition: Mesh.cpp:706
Definition: ColorToLabelMap.h:27
Mesh & loadFromFile(std::string const &meshFilename, std::string const &labelFilename)
Load mesh and labels from ASCII file using format auto detection.
Definition: Mesh.cpp:100
This finle contains definitions for label to color map types.
std::function< Label(Scalar, Scalar, Scalar)> ColorToLabelMap
Color to label map.
Definition: ColorToLabelMap.h:59
unsigned int Label
Label type.
Definition: Mesh.h:51
A basic ray datatype.
Definition: Ray.h:22
float Scalar
Scalar type of the vector space the embedding.
Definition: Mesh.h:42
Size triangleCount() const noexcept
Number of triangles.
Definition: Mesh.h:89
Definition: matlab-addons.c:7
This header contains the definition of the BarycentricPoint data type.
Definition: RayMeshIntersection.h:21
std::vector< std::array< T, 3 > > cartesianRepresentation(std::vector< BarycentricPoint< T, Index >> const &points) const
Converts a sequence of barycentric coordinates into caresian coordinates.
Definition: Mesh.h:224
Size vertexCount() const noexcept
Number of vertices.
Definition: Mesh.h:86
This file contains the definition of the RayMeshIntersection type.
auto withUnsafeVertexNormalPointer(F &&f) const noexcept(noexcept(f(std::declval< const NormalData >().data())))
Calls the given functional with an unsafe pointer to the raw per vertex normal storage.
Definition: Mesh.h:420
auto withUnsafeLabelPointer(F &&f) const noexcept(noexcept(f(LabelData().data())))
Calls the given functional with an unsafe pointer to the raw label storage.
Definition: Mesh.h:395
RayMeshIntersection< T > rayIntersection(Ray< T > const &ray) const
Computes intersection with the given ray and the mesh.
Definition: Mesh.cpp:550
Mesh< T > & upsample(std::size_t nTimes=1)
Upsample the mesh without touching the original vertices.
Definition: Mesh.cpp:627
Mesh(std::size_t nVertices, std::size_t nTriangles) noexcept(noexcept(VertexData(3 *nVertices, T{0}))&&noexcept(IndexData(3 *nTriangles, Index{0}))&&noexcept(LabelData(nVertices, Label{0}))&&noexcept(NormalData(3 *nVertices, T{0})))
Constructs an uninitialized mesh with the given vertex and triangle count.
Definition: Mesh.h:72
Definition of Ray data structure.
Mesh() noexcept(noexcept(VertexData())&&noexcept(IndexData())&&noexcept(LabelData())&&noexcept(VertexData()))
Constructs an empty mesh.
Definition: Mesh.h:68
BarycentricPoint data type. Represents a points on triangulation in barycentric coordinates.
Definition: BarycentricPoint.h:26
bool isEmpty() const noexcept
true iff the mesh is empty
Definition: Mesh.h:94