libigl v2.5.0
Loading...
Searching...
No Matches
CSGTree.h
Go to the documentation of this file.
1// This file is part of libigl, a simple c++ geometry processing library.
2//
3// Copyright (C) 2015 Alec Jacobson <alecjacobson@gmail.com>
4//
5// This Source Code Form is subject to the terms of the Mozilla Public License
6// v. 2.0. If a copy of the MPL was not distributed with this file, You can
7// obtain one at http://mozilla.org/MPL/2.0/.
8#ifndef IGL_COPYLEFT_CGAL_CSG_TREE_H
9#define IGL_COPYLEFT_CGAL_CSG_TREE_H
10
11#include "../../MeshBooleanType.h"
13#include "mesh_boolean.h"
14#include <CGAL/Exact_predicates_exact_constructions_kernel.h>
15#include <CGAL/number_utils.h>
16
17namespace igl
18{
19 namespace copyleft
20 {
21 namespace cgal
22 {
26 class CSGTree
27 {
28 public:
29 typedef CGAL::Epeck::FT ExactScalar;
30 //typedef Eigen::PlainObjectBase<DerivedF> POBF;
31 typedef Eigen::MatrixXi POBF;
32 typedef Eigen::Matrix<ExactScalar,Eigen::Dynamic,3> MatrixX3E;
33 typedef Eigen::VectorXi VectorJ;
34 private:
36 MatrixX3E m_V;
38 POBF m_F;
40 VectorJ m_J;
43 size_t m_number_of_birth_faces;
44 public:
46 {
47 }
48 //typedef Eigen::MatrixXd MatrixX3E;
49 //typedef Eigen::MatrixXi POBF;
50 // http://stackoverflow.com/a/3279550/148668
51 CSGTree(const CSGTree & other)
52 :
53 // copy things
54 m_V(other.m_V),
55 // This is an issue if m_F is templated
56 // https://forum.kde.org/viewtopic.php?f=74&t=128414
57 m_F(other.m_F),
58 m_J(other.m_J),
59 m_number_of_birth_faces(other.m_number_of_birth_faces)
60 {
61 }
62 // copy-swap idiom
63 friend void swap(CSGTree& first, CSGTree& second)
64 {
65 using std::swap;
66 // swap things
67 swap(first.m_V,second.m_V);
68 // This is an issue if m_F is templated, similar to
69 // https://forum.kde.org/viewtopic.php?f=74&t=128414
70 swap(first.m_F,second.m_F);
71 swap(first.m_J,second.m_J);
72 swap(first.m_number_of_birth_faces,second.m_number_of_birth_faces);
73 }
74 // Pass-by-value (aka copy)
76 {
77 swap(*this,other);
78 return *this;
79 }
80 CSGTree(CSGTree&& other):
81 // initialize via default constructor
82 CSGTree()
83 {
84 swap(*this,other);
85 }
92 const CSGTree & A,
93 const CSGTree & B,
94 const MeshBooleanType & type)
95 {
96 // conduct boolean operation
97 mesh_boolean(A.V(),A.F(),B.V(),B.F(),type,m_V,m_F,m_J);
98 // reindex m_J
99 std::for_each(m_J.data(),m_J.data()+m_J.size(),
100 [&](typename VectorJ::Scalar & j) -> void
101 {
102 if(j < A.F().rows())
103 {
104 j = A.J()(j);
105 }else
106 {
107 assert(j<(A.F().rows()+B.F().rows()));
108 j = A.number_of_birth_faces()+(B.J()(j-A.F().rows()));
109 }
110 });
111 m_number_of_birth_faces =
112 A.number_of_birth_faces() + B.number_of_birth_faces();
113 }
116 const CSGTree & A,
117 const CSGTree & B,
118 const std::string & s):
120 {
121 // do nothing (all done in constructor).
122 }
128 template <typename DerivedV>
129 CSGTree(const Eigen::PlainObjectBase<DerivedV> & V, const POBF & F)//:
130 // Possible Eigen bug:
131 // https://forum.kde.org/viewtopic.php?f=74&t=128414
132 //m_V(V.template cast<ExactScalar>()),m_F(F)
133 {
134 m_V = V.template cast<ExactScalar>();
135 m_F = F;
136 // number of faces
137 m_number_of_birth_faces = m_F.rows();
138 // identity birth index
139 m_J = VectorJ::LinSpaced(
140 m_number_of_birth_faces,0,m_number_of_birth_faces-1);
141 }
142 // Returns reference to resulting mesh vertices m_V in exact scalar
143 // representation
144 const MatrixX3E & V() const
145 {
146 return m_V;
147 }
148 // Returns mesh vertices in the desired output type, casting when
149 // appropriate to floating precision.
150 template <typename DerivedV>
151 DerivedV cast_V() const
152 {
153 DerivedV dV;
154 dV.resize(m_V.rows(),m_V.cols());
155 for(int i = 0;i<m_V.rows();i++)
156 {
157 for(int j = 0;j<m_V.cols();j++)
158 {
159 dV(i,j) = CGAL::to_double(m_V(i,j));
160 }
161 }
162 return dV;
163 }
164 // Returns reference to resulting mesh faces m_F
165 const POBF & F() const
166 {
167 return m_F;
168 }
169 // Returns reference to "birth parents" indices into [F1;F2;...;Fn]
170 // where F1, ... , Fn are the face lists of the leaf ("original") input
171 // meshes.
172 const VectorJ & J() const
173 {
174 return m_J;
175 }
176 // The number of leaf faces = #F1 + #F2 + ... + #Fn
177 const size_t & number_of_birth_faces() const
178 {
179 return m_number_of_birth_faces;
180 }
181 };
182 }
183 }
184}
185
186
187#endif
Class for defining and computing a constructive solid geometry result out of a tree of boolean operat...
Definition CSGTree.h:27
const size_t & number_of_birth_faces() const
Definition CSGTree.h:177
CSGTree(const CSGTree &A, const CSGTree &B, const std::string &s)
This is an overloaded member function, provided for convenience. It differs from the above function o...
Definition CSGTree.h:115
CGAL::Epeck::FT ExactScalar
Definition CSGTree.h:29
CSGTree(const Eigen::PlainObjectBase< DerivedV > &V, const POBF &F)
"Leaf" node with identity operation on assumed "solid" mesh (V,F)
Definition CSGTree.h:129
CSGTree & operator=(CSGTree other)
Definition CSGTree.h:75
const POBF & F() const
Definition CSGTree.h:165
friend void swap(CSGTree &first, CSGTree &second)
Definition CSGTree.h:63
CSGTree(CSGTree &&other)
Definition CSGTree.h:80
CSGTree(const CSGTree &other)
Definition CSGTree.h:51
DerivedV cast_V() const
Definition CSGTree.h:151
const MatrixX3E & V() const
Definition CSGTree.h:144
Eigen::VectorXi VectorJ
Definition CSGTree.h:33
CSGTree()
Definition CSGTree.h:45
CSGTree(const CSGTree &A, const CSGTree &B, const MeshBooleanType &type)
Construct and compute a boolean operation on existing CSGTree nodes.
Definition CSGTree.h:91
Eigen::MatrixXi POBF
Definition CSGTree.h:31
const VectorJ & J() const
Definition CSGTree.h:172
Eigen::Matrix< ExactScalar, Eigen::Dynamic, 3 > MatrixX3E
Definition CSGTree.h:32
bool mesh_boolean(const Eigen::MatrixBase< DerivedVA > &VA, const Eigen::MatrixBase< DerivedFA > &FA, const Eigen::MatrixBase< DerivedVB > &VB, const Eigen::MatrixBase< DerivedFB > &FB, const MeshBooleanType &type, Eigen::PlainObjectBase< DerivedVC > &VC, Eigen::PlainObjectBase< DerivedFC > &FC, Eigen::PlainObjectBase< DerivedJ > &J)
Compute Boolean csg operations on "solid", consistently oriented meshes.
bool string_to_mesh_boolean_type(const std::string &s, MeshBooleanType &type)
Convert string to boolean type.
Definition AABB.h:17
MeshBooleanType
Boolean operation types.
Definition MeshBooleanType.h:14