preCICE v3.1.1
Loading...
Searching...
No Matches
Tetrahedron.cpp
Go to the documentation of this file.
1#include "Tetrahedron.hpp"
2#include <Eigen/Core>
3#include <algorithm>
4#include <boost/concept/assert.hpp>
6#include "math/geometry.hpp"
7#include "mesh/Vertex.hpp"
8#include "utils/EigenIO.hpp"
9
10namespace precice::mesh {
11
13 Vertex &vertexOne,
14 Vertex &vertexTwo,
15 Vertex &vertexThree,
16 Vertex &vertexFour)
17 : _vertices({&vertexOne, &vertexTwo, &vertexThree, &vertexFour})
18{
19 PRECICE_ASSERT(vertexOne.getDimensions() == vertexTwo.getDimensions(),
20 vertexOne.getDimensions(), vertexTwo.getDimensions());
21 PRECICE_ASSERT(vertexOne.getDimensions() == vertexThree.getDimensions(),
22 vertexOne.getDimensions(), vertexThree.getDimensions());
23 PRECICE_ASSERT(vertexOne.getDimensions() == vertexFour.getDimensions(),
24 vertexOne.getDimensions(), vertexFour.getDimensions());
25 PRECICE_ASSERT(getDimensions() == 3, getDimensions());
26
28 (&vertexOne != &vertexTwo) &&
29 (&vertexOne != &vertexThree) &&
30 (&vertexOne != &vertexFour) &&
31 (&vertexTwo != &vertexThree) &&
32 (&vertexTwo != &vertexFour) &&
33 (&vertexThree != &vertexFour),
34 "Tetrahedron vertices are not unique!");
35
36 std::sort(_vertices.begin(), _vertices.end(),
37 [](const Vertex *lhs, const Vertex *rhs) { return *lhs < *rhs; });
38}
39
41{
42 return math::geometry::tetraVolume(vertex(0).getCoords(), vertex(1).getCoords(), vertex(2).getCoords(), vertex(3).getCoords());
43}
44
46{
47 return _vertices[0]->getDimensions();
48}
49
50const Eigen::VectorXd Tetrahedron::getCenter() const
51{
52 return (vertex(0).getCoords() + vertex(1).getCoords() + vertex(2).getCoords() + vertex(3).getCoords()) / 4.0;
53}
54
56{
57 auto center = getCenter();
58 return std::max({(center - vertex(0).getCoords()).norm(),
59 (center - vertex(1).getCoords()).norm(),
60 (center - vertex(2).getCoords()).norm(),
61 (center - vertex(3).getCoords()).norm()});
62}
63
64bool Tetrahedron::operator==(const Tetrahedron &other) const
65{
66 return std::is_permutation(_vertices.begin(), _vertices.end(), other._vertices.begin(),
67 [](const Vertex *v1, const Vertex *v2) { return *v1 == *v2; });
68}
69
70bool Tetrahedron::operator!=(const Tetrahedron &other) const
71{
72 return !(*this == other);
73}
74
76{
77 // Show 6 edges: 0-1, 0-2, 0-3, 1-2, 1-3, 2_3
79 const auto &v0 = t.vertex(0);
80 const auto &v1 = t.vertex(1);
81 const auto &v2 = t.vertex(2);
82 const auto &v3 = t.vertex(3);
83
84 return os << "MULTILINESTRING ("
85 << "(" << v0.getCoords().transpose().format(wkt()) << ", " << v1.getCoords().transpose().format(wkt()) << "), "
86 << "(" << v0.getCoords().transpose().format(wkt()) << ", " << v2.getCoords().transpose().format(wkt()) << "), "
87 << "(" << v0.getCoords().transpose().format(wkt()) << ", " << v3.getCoords().transpose().format(wkt()) << "), "
88 << "(" << v1.getCoords().transpose().format(wkt()) << ", " << v2.getCoords().transpose().format(wkt()) << "), "
89 << "(" << v1.getCoords().transpose().format(wkt()) << ", " << v3.getCoords().transpose().format(wkt()) << "), "
90 << "(" << v2.getCoords().transpose().format(wkt()) << ", " << v3.getCoords().transpose().format(wkt()) << ")"
91 << ")";
92}
93
94} // namespace precice::mesh
#define PRECICE_ASSERT(...)
Definition assertion.hpp:87
Tetrahedron of a mesh, defined by 4 vertices.
Vertex & vertex(int i)
Returns tetrahedron vertex with index 0, 1, 2 or 3.
std::array< Vertex *, 4 > _vertices
Vertices defining the Tetrahedron.
const Eigen::VectorXd getCenter() const
Returns the barycenter of the tetrahedron.
double getEnclosingRadius() const
Returns the radius of the sphere enclosing the tetrahedron.
double getVolume() const
Returns the unsigned volume of the tetrahedron.
bool operator!=(const Tetrahedron &other) const
Not equal, implemented in terms of equal.
bool operator==(const Tetrahedron &other) const
Compares two Tetrahedrons for equality.
int getDimensions() const
Returns dimensionalty of space the Tetrahedron is embedded in.
Tetrahedron(Vertex &vertexOne, Vertex &vertexTwo, Vertex &vertexThree, Vertex &vertexFour)
Vertex of a mesh.
Definition Vertex.hpp:16
Eigen::VectorXd getCoords() const
Returns the coordinates of the vertex.
Definition Vertex.hpp:116
T is_permutation(T... args)
T max(T... args)
double tetraVolume(const Eigen::Vector3d &a, const Eigen::Vector3d &b, const Eigen::Vector3d &c, const Eigen::Vector3d &d)
Computes the (unsigned) area of a triangle in 3D.
Definition geometry.cpp:116
provides Mesh, Data and primitives.
std::ostream & operator<<(std::ostream &os, const BoundingBox &bb)
Eigen::IOFormat wkt()
Definition EigenIO.hpp:9
T sort(T... args)