PolyFEM
Loading...
Searching...
No Matches
Mesh.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <polyfem/Common.hpp>
8
9#include <Eigen/Dense>
10#include <geogram/mesh/mesh.h>
11
12#include <memory>
13
14namespace polyfem
15{
16 namespace mesh
17 {
36
38 class Mesh
39 {
40 protected:
43 {
44 public:
45 int v1, v2;
46 Eigen::MatrixXd nodes;
47 };
48
51 {
52 public:
53 int v1, v2, v3;
54 Eigen::MatrixXd nodes;
55 };
56
59 {
60 public:
61 int v1, v2, v3, v4;
62 Eigen::MatrixXd nodes;
63 };
64
65 public:
72 static std::unique_ptr<Mesh> create(const std::string &path, const bool non_conforming = false);
73
80 static std::unique_ptr<Mesh> create(GEO::Mesh &M, const bool non_conforming = false);
81
89 static std::unique_ptr<Mesh> create(const Eigen::MatrixXd &vertices, const Eigen::MatrixXi &cells, const bool non_conforming = false);
90
97 static std::unique_ptr<Mesh> create(const int dim, const bool non_conforming = false);
98
101 virtual std::unique_ptr<Mesh> copy() const = 0;
102
103 protected:
107 Mesh() = default;
108
109 public:
112 virtual ~Mesh() = default;
116 Mesh(Mesh &&) = default;
121 Mesh &operator=(Mesh &&) = default;
125 Mesh(const Mesh &) = default;
130 Mesh &operator=(const Mesh &) = default;
131
137 virtual void refine(const int n_refinement, const double t) = 0;
138
143 virtual bool is_volume() const = 0;
148 int dimension() const { return (is_volume() ? 3 : 2); }
153 virtual bool is_conforming() const = 0;
158 int n_elements() const { return (is_volume() ? n_cells() : n_faces()); }
163 int n_boundary_elements() const { return (is_volume() ? n_faces() : n_edges()); }
164
169 virtual int n_cells() const = 0;
174 virtual int n_faces() const = 0;
178 virtual int n_edges() const = 0;
182 virtual int n_vertices() const = 0;
183
188 virtual int n_face_vertices(const int f_id) const = 0;
193 virtual int n_cell_vertices(const int c_id) const = 0;
199 virtual int edge_vertex(const int e_id, const int lv_id) const = 0;
205 virtual int face_vertex(const int f_id, const int lv_id) const = 0;
211 virtual int cell_vertex(const int f_id, const int lv_id) const = 0;
217 int element_vertex(const int el_id, const int lv_id) const
218 {
219 return (is_volume() ? cell_vertex(el_id, lv_id) : face_vertex(el_id, lv_id));
220 }
221
222 int boundary_element_vertex(const int primitive_id, const int lv_id) const
223 {
224 return (is_volume() ? face_vertex(primitive_id, lv_id) : edge_vertex(primitive_id, lv_id));
225 }
226
231 virtual bool is_boundary_vertex(const int vertex_global_id) const = 0;
236 virtual bool is_boundary_edge(const int edge_global_id) const = 0;
241 virtual bool is_boundary_face(const int face_global_id) const = 0;
246 virtual bool is_boundary_element(const int element_global_id) const = 0;
247
248 virtual bool save(const std::string &path) const = 0;
249
250 private:
256 virtual bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F) = 0;
257
258 public:
263 virtual void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector<std::vector<int>> &nodes) = 0;
267 inline const Eigen::MatrixXi &orders() const { return orders_; }
271 inline bool is_rational() const { return is_rational_; }
275 inline void set_is_rational(const bool in_is_rational) { is_rational_ = in_is_rational; }
276
279 virtual void normalize() = 0;
280
282 virtual void compute_elements_tag() = 0;
284 virtual void update_elements_tag() { assert(false); }
285
290 virtual double edge_length(const int gid) const
291 {
292 assert(false);
293 return 0;
294 }
299 virtual double quad_area(const int gid) const
300 {
301 assert(false);
302 return 0;
303 }
308 virtual double tri_area(const int gid) const
309 {
310 assert(false);
311 return 0;
312 }
313
318 virtual RowVectorNd point(const int global_index) const = 0;
323 virtual void set_point(const int global_index, const RowVectorNd &p) = 0;
324
329 virtual RowVectorNd edge_barycenter(const int e) const = 0;
334 virtual RowVectorNd face_barycenter(const int f) const = 0;
339 virtual RowVectorNd cell_barycenter(const int c) const = 0;
340
344 void edge_barycenters(Eigen::MatrixXd &barycenters) const;
348 void face_barycenters(Eigen::MatrixXd &barycenters) const;
352 void cell_barycenters(Eigen::MatrixXd &barycenters) const;
356 virtual void compute_element_barycenters(Eigen::MatrixXd &barycenters) const = 0;
357
361 virtual void elements_boxes(std::vector<std::array<Eigen::Vector3d, 2>> &boxes) const = 0;
367 virtual void barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const = 0;
368
373 virtual void bounding_box(RowVectorNd &min, RowVectorNd &max) const = 0;
374
379 bool is_spline_compatible(const int el_id) const;
384 bool is_cube(const int el_id) const;
389 bool is_polytope(const int el_id) const;
394 bool is_simplex(const int el_id) const;
395
399 const std::vector<ElementType> &elements_tag() const { return elements_tag_; }
404 void set_tag(const int el, const ElementType type) { elements_tag_[el] = type; }
405
409 void compute_node_ids(const std::function<int(const size_t, const RowVectorNd &, bool)> &marker);
410
414 virtual void load_boundary_ids(const std::string &path);
419 virtual void compute_boundary_ids(const double eps) = 0;
423 virtual void compute_boundary_ids(const std::function<int(const RowVectorNd &)> &marker) = 0;
427 virtual void compute_boundary_ids(const std::function<int(const RowVectorNd &, bool)> &marker) = 0;
431 virtual void compute_boundary_ids(const std::function<int(const size_t, const RowVectorNd &, bool)> &marker) = 0;
435 virtual void compute_boundary_ids(const std::function<int(const std::vector<int> &, bool)> &marker) = 0;
439 virtual void compute_boundary_ids(const std::function<int(const size_t, const std::vector<int> &, const RowVectorNd &, bool)> &marker) = 0;
440
444 virtual void compute_body_ids(const std::function<int(const size_t, const RowVectorNd &)> &marker) = 0;
448 virtual void set_boundary_ids(const std::vector<int> &boundary_ids)
449 {
450 assert(boundary_ids.size() == n_boundary_elements());
451 boundary_ids_ = boundary_ids;
452 }
456 virtual void set_body_ids(const std::vector<int> &body_ids)
457 {
458 assert(body_ids.size() == n_elements());
459 body_ids_ = body_ids;
460 }
461
466 virtual int get_default_boundary_id(const int primitive) const
467 {
468 if (is_volume() ? is_boundary_face(primitive) : is_boundary_edge(primitive))
469 return std::numeric_limits<int>::max(); // default for no selected boundary
470 else
471 return -1; // default for no boundary
472 }
473
478 virtual int get_boundary_id(const int primitive) const
479 {
480 return has_boundary_ids() ? (boundary_ids_.at(primitive)) : get_default_boundary_id(primitive);
481 }
482
487 virtual int get_node_id(const int node_id) const
488 {
489 if (has_node_ids())
490 return node_ids_.at(node_id);
491 else
492 return -1; // default for no boundary
493 }
494
498 void update_nodes(const Eigen::VectorXi &in_node_to_node);
499
504 virtual int get_body_id(const int primitive) const
505 {
506 if (has_body_ids())
507 return body_ids_.at(primitive);
508 else
509 return 0;
510 }
513 virtual const std::vector<int> &get_body_ids() const
514 {
515 return body_ids_;
516 }
520 bool has_node_ids() const { return !node_ids_.empty(); }
524 bool has_boundary_ids() const { return !boundary_ids_.empty(); }
528 virtual bool has_body_ids() const { return !body_ids_.empty(); }
529
534 virtual void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1) const = 0;
540 virtual void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1, const std::vector<bool> &valid_elements) const = 0;
546 virtual void triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector<int> &ranges) const = 0;
547
552 const std::vector<double> &cell_weights(const int cell_index) const { return cell_weights_[cell_index]; }
556 void set_cell_weights(const std::vector<std::vector<double>> &in_cell_weights) { cell_weights_ = in_cell_weights; }
557
560 virtual void prepare_mesh(){};
561
565 bool has_poly() const
566 {
567 for (int i = 0; i < n_elements(); ++i)
568 {
569 if (is_polytope(i))
570 return true;
571 }
572
573 return false;
574 }
575
579 bool is_simplicial() const
580 {
581 for (int i = 0; i < n_elements(); ++i)
582 {
583 if (!is_simplex(i))
584 return false;
585 }
586
587 return true;
588 }
589
593 inline bool is_linear() const { return orders_.size() == 0 || orders_.maxCoeff() == 1; }
594
598 std::vector<std::pair<int, int>> edges() const;
602 std::vector<std::vector<int>> faces() const;
603
607 std::unordered_map<std::pair<int, int>, size_t, polyfem::utils::HashPair> edges_to_ids() const;
611 std::unordered_map<std::vector<int>, size_t, polyfem::utils::HashVector> faces_to_ids() const;
612
616 inline const Eigen::VectorXi &in_ordered_vertices() const { return in_ordered_vertices_; }
620 inline const Eigen::MatrixXi &in_ordered_edges() const { return in_ordered_edges_; }
624 inline const Eigen::MatrixXi &in_ordered_faces() const { return in_ordered_faces_; }
625
629 virtual void append(const Mesh &mesh);
630
634 void append(const std::unique_ptr<Mesh> &mesh)
635 {
636 if (mesh != nullptr)
637 append(*mesh);
638 }
639
643 void apply_affine_transformation(const MatrixNd &A, const VectorNd &b);
644
645 protected:
650 virtual bool load(const std::string &path) = 0;
655 virtual bool load(const GEO::Mesh &M) = 0;
656
658 std::vector<ElementType> elements_tag_;
660 std::vector<int> node_ids_;
662 std::vector<int> boundary_ids_;
664 std::vector<int> body_ids_;
666 Eigen::MatrixXi orders_;
668 bool is_rational_ = false;
669
671 std::vector<EdgeNodes> edge_nodes_;
673 std::vector<FaceNodes> face_nodes_;
675 std::vector<CellNodes> cell_nodes_;
677 std::vector<std::vector<double>> cell_weights_;
678
680 Eigen::VectorXi in_ordered_vertices_;
682 Eigen::MatrixXi in_ordered_edges_;
684 Eigen::MatrixXi in_ordered_faces_;
685 };
686 } // namespace mesh
687} // namespace polyfem
int V
Class to store the high-order cells nodes.
Definition Mesh.hpp:59
Eigen::MatrixXd nodes
Definition Mesh.hpp:62
Class to store the high-order edge nodes.
Definition Mesh.hpp:43
Eigen::MatrixXd nodes
Definition Mesh.hpp:46
Class to store the high-order face nodes.
Definition Mesh.hpp:51
Eigen::MatrixXd nodes
Definition Mesh.hpp:54
Abstract mesh class to capture 2d/3d conforming and non-conforming meshes.
Definition Mesh.hpp:39
int n_elements() const
utitlity to return the number of elements, cells or faces in 3d and 2d
Definition Mesh.hpp:158
virtual int n_vertices() const =0
number of vertices
virtual int get_body_id(const int primitive) const
Get the volume selection of an element (cell in 3d, face in 2d)
Definition Mesh.hpp:504
virtual void compute_element_barycenters(Eigen::MatrixXd &barycenters) const =0
utility for 2d/3d.
virtual double edge_length(const int gid) const
edge length
Definition Mesh.hpp:290
void set_cell_weights(const std::vector< std::vector< double > > &in_cell_weights)
Set the cell weights for rational polynomial meshes.
Definition Mesh.hpp:556
const Eigen::MatrixXi & in_ordered_edges() const
Order of the input edges.
Definition Mesh.hpp:620
bool is_polytope(const int el_id) const
checks if element is polygon compatible
Definition Mesh.cpp:361
Eigen::MatrixXi orders_
list of geometry orders, one per cell
Definition Mesh.hpp:666
virtual void compute_boundary_ids(const std::function< int(const size_t, const RowVectorNd &, bool)> &marker)=0
computes boundary selections based on a function
std::unordered_map< std::pair< int, int >, size_t, polyfem::utils::HashPair > edges_to_ids() const
map from edge (pair of v id) to the id of the edge
Definition Mesh.cpp:455
bool is_rational() const
check if curved mesh has rational polynomials elements
Definition Mesh.hpp:271
Mesh()=default
Construct a new Mesh object.
std::vector< ElementType > elements_tag_
list of element types
Definition Mesh.hpp:658
virtual RowVectorNd edge_barycenter(const int e) const =0
edge barycenter
virtual void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1) const =0
Get all the edges.
int element_vertex(const int el_id, const int lv_id) const
id of the vertex of a element
Definition Mesh.hpp:217
bool is_rational_
stores if the mesh is rational
Definition Mesh.hpp:668
bool has_boundary_ids() const
checks if surface selections are available
Definition Mesh.hpp:524
virtual void compute_boundary_ids(const std::function< int(const size_t, const std::vector< int > &, const RowVectorNd &, bool)> &marker)=0
computes boundary selections based on a function
virtual double tri_area(const int gid) const
area of a tri face of a tet mesh
Definition Mesh.hpp:308
void cell_barycenters(Eigen::MatrixXd &barycenters) const
all cells barycenters
Definition Mesh.cpp:316
virtual void get_edges(Eigen::MatrixXd &p0, Eigen::MatrixXd &p1, const std::vector< bool > &valid_elements) const =0
Get all the edges according to valid_elements selection.
virtual void compute_boundary_ids(const std::function< int(const std::vector< int > &, bool)> &marker)=0
computes boundary selections based on a function
bool is_simplicial() const
checks if the mesh is simplicial
Definition Mesh.hpp:579
virtual bool is_conforming() const =0
if the mesh is conforming
virtual void bounding_box(RowVectorNd &min, RowVectorNd &max) const =0
computes the bbox of the mesh
virtual void compute_boundary_ids(const std::function< int(const RowVectorNd &, bool)> &marker)=0
computes boundary selections based on a function
virtual void barycentric_coords(const RowVectorNd &p, const int el_id, Eigen::MatrixXd &coord) const =0
constructs barycentric coodiantes for a point p.
virtual RowVectorNd face_barycenter(const int f) const =0
face barycenter
virtual void update_elements_tag()
Update elements types.
Definition Mesh.hpp:284
Mesh & operator=(const Mesh &)=default
virtual void compute_elements_tag()=0
compute element types, see ElementType
virtual bool build_from_matrices(const Eigen::MatrixXd &V, const Eigen::MatrixXi &F)=0
build a mesh from matrices
virtual void set_point(const int global_index, const RowVectorNd &p)=0
Set the point.
bool is_cube(const int el_id) const
checks if element is cube compatible
Definition Mesh.cpp:348
bool has_node_ids() const
checks if points selections are available
Definition Mesh.hpp:520
virtual void set_body_ids(const std::vector< int > &body_ids)
Set the volume sections.
Definition Mesh.hpp:456
virtual RowVectorNd point(const int global_index) const =0
point coordinates
const Eigen::MatrixXi & orders() const
order of each element
Definition Mesh.hpp:267
virtual bool is_boundary_face(const int face_global_id) const =0
is face boundary
Eigen::MatrixXi in_ordered_faces_
Order of the input faces, TODO: change to std::vector of Eigen::Vector.
Definition Mesh.hpp:684
void compute_node_ids(const std::function< int(const size_t, const RowVectorNd &, bool)> &marker)
computes boundary selections based on a function
Definition Mesh.cpp:383
virtual int get_boundary_id(const int primitive) const
Get the boundary selection of an element (face in 3d, edge in 2d)
Definition Mesh.hpp:478
virtual void triangulate_faces(Eigen::MatrixXi &tris, Eigen::MatrixXd &pts, std::vector< int > &ranges) const =0
generate a triangular representation of every face
virtual bool is_boundary_vertex(const int vertex_global_id) const =0
is vertex boundary
bool is_simplex(const int el_id) const
checks if element is simples compatible
Definition Mesh.cpp:418
void face_barycenters(Eigen::MatrixXd &barycenters) const
all face barycenters
Definition Mesh.cpp:307
virtual double quad_area(const int gid) const
area of a quad face of an hex mesh
Definition Mesh.hpp:299
virtual bool has_body_ids() const
checks if volumes selections are available
Definition Mesh.hpp:528
virtual void normalize()=0
normalize the mesh
bool is_spline_compatible(const int el_id) const
checks if element is spline compatible
Definition Mesh.cpp:328
virtual void prepare_mesh()
method used to finalize the mesh.
Definition Mesh.hpp:560
virtual void load_boundary_ids(const std::string &path)
loads the boundary selections for a file
Definition Mesh.cpp:395
void set_is_rational(const bool in_is_rational)
Set the is rational object.
Definition Mesh.hpp:275
Mesh & operator=(Mesh &&)=default
Copy constructor.
std::vector< int > boundary_ids_
list of surface labels
Definition Mesh.hpp:662
bool is_linear() const
check if the mesh is linear
Definition Mesh.hpp:593
void apply_affine_transformation(const MatrixNd &A, const VectorNd &b)
Apply an affine transformation to the vertex positions .
Definition Mesh.cpp:628
std::vector< int > node_ids_
list of node labels
Definition Mesh.hpp:660
std::unordered_map< std::vector< int >, size_t, polyfem::utils::HashVector > faces_to_ids() const
map from face (tuple of v id) to the id of the face
Definition Mesh.cpp:471
virtual void compute_boundary_ids(const double eps)=0
computes the selection based on the bbx of the mesh.
virtual void refine(const int n_refinement, const double t)=0
refine the mesh
std::vector< CellNodes > cell_nodes_
high-order nodes associates to cells
Definition Mesh.hpp:675
std::vector< std::vector< double > > cell_weights_
weights associates to cells for rational polynomail meshes
Definition Mesh.hpp:677
const std::vector< double > & cell_weights(const int cell_index) const
weights for rational polynomial meshes
Definition Mesh.hpp:552
virtual bool load(const GEO::Mesh &M)=0
loads a mesh from a geo mesh
virtual bool is_volume() const =0
checks if mesh is volume
std::vector< std::pair< int, int > > edges() const
list of sorted edges.
Definition Mesh.cpp:423
void update_nodes(const Eigen::VectorXi &in_node_to_node)
Update the node ids to reorder them.
Definition Mesh.cpp:367
virtual int edge_vertex(const int e_id, const int lv_id) const =0
id of the edge vertex
virtual std::unique_ptr< Mesh > copy() const =0
Create a copy of the mesh.
virtual void compute_boundary_ids(const std::function< int(const RowVectorNd &)> &marker)=0
computes boundary selections based on a function
void append(const std::unique_ptr< Mesh > &mesh)
appends a new mesh to the end of this, utility that takes pointer, calls other one
Definition Mesh.hpp:634
virtual bool is_boundary_edge(const int edge_global_id) const =0
is edge boundary
bool has_poly() const
checks if the mesh has polytopes
Definition Mesh.hpp:565
const Eigen::VectorXi & in_ordered_vertices() const
Order of the input vertices.
Definition Mesh.hpp:616
std::vector< int > body_ids_
list of volume labels
Definition Mesh.hpp:664
static std::unique_ptr< Mesh > create(const std::string &path, const bool non_conforming=false)
factory to build the proper mesh
Definition Mesh.cpp:173
virtual int get_default_boundary_id(const int primitive) const
Get the default boundary selection of an element (face in 3d, edge in 2d)
Definition Mesh.hpp:466
int dimension() const
utily for dimension
Definition Mesh.hpp:148
virtual int n_cells() const =0
number of cells
virtual const std::vector< int > & get_body_ids() const
Get the volume selection of all elements (cells in 3d, faces in 2d)
Definition Mesh.hpp:513
std::vector< FaceNodes > face_nodes_
high-order nodes associates to faces
Definition Mesh.hpp:673
virtual bool load(const std::string &path)=0
loads a mesh from the path
virtual int n_faces() const =0
number of faces
const std::vector< ElementType > & elements_tag() const
Returns the elements types.
Definition Mesh.hpp:399
std::vector< EdgeNodes > edge_nodes_
high-order nodes associates to edges
Definition Mesh.hpp:671
std::vector< std::vector< int > > faces() const
list of sorted faces.
Definition Mesh.cpp:439
virtual void attach_higher_order_nodes(const Eigen::MatrixXd &V, const std::vector< std::vector< int > > &nodes)=0
attach high order nodes
int n_boundary_elements() const
utitlity to return the number of boundary elements, faces or edges in 3d and 2d
Definition Mesh.hpp:163
Eigen::MatrixXi in_ordered_edges_
Order of the input edges.
Definition Mesh.hpp:682
void set_tag(const int el, const ElementType type)
changes the element type
Definition Mesh.hpp:404
virtual void append(const Mesh &mesh)
appends a new mesh to the end of this
Definition Mesh.cpp:490
virtual void compute_body_ids(const std::function< int(const size_t, const RowVectorNd &)> &marker)=0
computes boundary selections based on a function
Mesh(Mesh &&)=default
Construct a new Mesh object.
virtual bool save(const std::string &path) const =0
virtual RowVectorNd cell_barycenter(const int c) const =0
cell barycenter
const Eigen::MatrixXi & in_ordered_faces() const
Order of the input edges.
Definition Mesh.hpp:624
virtual int n_edges() const =0
number of edges
virtual void set_boundary_ids(const std::vector< int > &boundary_ids)
Set the boundary selection from a vector.
Definition Mesh.hpp:448
void edge_barycenters(Eigen::MatrixXd &barycenters) const
all edges barycenters
Definition Mesh.cpp:298
virtual int cell_vertex(const int f_id, const int lv_id) const =0
id of the vertex of a cell
virtual int n_face_vertices(const int f_id) const =0
number of vertices of a face
virtual void elements_boxes(std::vector< std::array< Eigen::Vector3d, 2 > > &boxes) const =0
constructs a box around every element (3d cell, 2d face)
Eigen::VectorXi in_ordered_vertices_
Order of the input vertices.
Definition Mesh.hpp:680
virtual bool is_boundary_element(const int element_global_id) const =0
is cell boundary
virtual int n_cell_vertices(const int c_id) const =0
number of vertices of a cell
virtual int get_node_id(const int node_id) const
Get the boundary selection of a node.
Definition Mesh.hpp:487
virtual ~Mesh()=default
Destroy the Mesh object.
virtual int face_vertex(const int f_id, const int lv_id) const =0
id of the face vertex
Mesh(const Mesh &)=default
Construct a new Mesh object.
int boundary_element_vertex(const int primitive_id, const int lv_id) const
Definition Mesh.hpp:222
ElementType
Type of Element, check [Poly-Spline Finite Element Method] for a complete description.
Definition Mesh.hpp:23
@ REGULAR_INTERIOR_CUBE
Triangle/tet element.
@ UNDEFINED
Boundary polytope.
@ REGULAR_BOUNDARY_CUBE
Quad/Hex incident to more than 1 singular vertices (should not happen in 2D)
@ MULTI_SINGULAR_BOUNDARY_CUBE
Quad incident to exactly 1 singular vertex (in 2D); hex incident to exactly 1 singular interior edge,...
@ MULTI_SINGULAR_INTERIOR_CUBE
Quad/hex incident to exactly 1 singular vertex (in 2D) or edge (in 3D)
@ SIMPLE_SINGULAR_INTERIOR_CUBE
Regular quad/hex inside a 3^n patch.
@ INTERFACE_CUBE
Boundary hex that is not regular nor SimpleSingularBoundaryCube.
@ INTERIOR_POLYTOPE
Quad/hex that is at the interface with a polytope (if a cube has both external boundary and and inter...
@ SIMPLE_SINGULAR_BOUNDARY_CUBE
Boundary quad/hex, where all boundary vertices/edges are incident to at most 2 quads/hexes.
@ BOUNDARY_POLYTOPE
Interior polytope.
Eigen::Matrix< double, Eigen::Dynamic, 1, 0, 3, 1 > VectorNd
Definition Types.hpp:9
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor, 3, 3 > MatrixNd
Definition Types.hpp:12
Eigen::Matrix< double, 1, Eigen::Dynamic, Eigen::RowMajor, 1, 3 > RowVectorNd
Definition Types.hpp:11