PeriDyno 1.2.1
Loading...
Searching...
No Matches
FastSweepingMethod.cpp
Go to the documentation of this file.
2
3namespace dyno
4{
6
7 template<typename TDataType>
10 {
11 this->varPadding()->setMin(2);
12 }
13
14 template<typename TDataType>
16 {
17 faceList.clear();
18 vertList.clear();
19 phi.clear();
20 }
21
22 template<typename TDataType>
24 {
25 this->loadClosedSurface();
26
27 if (this->outLevelSet()->isEmpty()) {
28 this->outLevelSet()->allocate();
29 }
30 DistanceField3D<TDataType>& sdf = this->outLevelSet()->getDataPtr()->getSDF();
31
32 Coord p0(origin.x, origin.y, origin.z);
33 Coord p1(maxPoint.x, maxPoint.y, maxPoint.z);
34
35 sdf.setSpace(p0, p1, this->varSpacing()->getValue());
36 sdf.setDistance(phi);
37 }
38
39 template<typename TDataType>
41 {
42 Vec3f min_box(std::numeric_limits<float>::max(), std::numeric_limits<float>::max(), std::numeric_limits<float>::max());
43 Vec3f max_box(-std::numeric_limits<float>::max(), -std::numeric_limits<float>::max(), -std::numeric_limits<float>::max());
44
45 auto& mPoints = inTriangleSet()->getDataPtr()->getPoints();
46 auto& mTriangles = inTriangleSet()->getDataPtr()->getTriangles();
47
48 CArray<Coord> cPoints;
49 cPoints.assign(mPoints);
50
52 cTriangles.resize(mTriangles.size());
53 cTriangles.assign(mTriangles);
54
55 vertList.clear();
56 faceList.clear();
57
58 for (int i = 0; i < cPoints.size(); i++) {
59 Coord point = cPoints[i];
60 vertList.pushBack(point);
61 min_box = min_box.minimum(point);
62 max_box = max_box.maximum(point);
63 }
64
65 for (int i = 0; i < cTriangles.size(); i++) {
66 faceList.pushBack(Vec3ui(cTriangles[i][0], cTriangles[i][1], cTriangles[i][2]));
67 }
68
69 Real dx = this->varSpacing()->getValue();
70 uint padding = this->varPadding()->getValue();
71
72 Vec3f unit(1, 1, 1);
73 min_box -= padding * dx * unit;
74 max_box += padding * dx * unit;
75
76 ni = std::ceil((max_box[0] - min_box[0]) / dx);
77 nj = std::ceil((max_box[1] - min_box[1]) / dx);
78 nk = std::ceil((max_box[2] - min_box[2]) / dx);
79
80 origin = min_box;
81 maxPoint = max_box;
82
84
85 cPoints.clear();
86 cTriangles.clear();
87 printf("Uniform grids: %f %f %f, %f %f %f, %f, %d %d %d, %d \n", origin[0], origin[1], origin[2], maxPoint[0], maxPoint[1], maxPoint[2], dx, ni, nj, nk, padding);
88 }
89
90 // find distance x0 is from segment x1-x2
91 static float point_segment_distance(const Vec3f &x0, const Vec3f &x1, const Vec3f &x2)
92 {
93 Vec3f dx(x2 - x1);
94 double m2 = dx.normSquared();
95 // find parameter value of closest point on segment
96 float s12 = (float)((x2 - x0).dot(dx) / m2);
97 if (s12 < 0) {
98 s12 = 0;
99 }
100 else if (s12 > 1) {
101 s12 = 1;
102 }
103 // and find the distance
104 return (x0 - (s12*x1 + (1 - s12)*x2)).norm();
105 }
106
107 // find distance x0 is from triangle x1-x2-x3
108 static float point_triangle_distance(const Vec3f &x0, const Vec3f &x1, const Vec3f &x2, const Vec3f &x3)
109 {
110 // first find barycentric coordinates of closest point on infinite plane
111 Vec3f x13(x1 - x3), x23(x2 - x3), x03(x0 - x3);
112 float m13 = x13.normSquared(), m23 = x23.normSquared(), d = x13.dot(x23);
113 float invdet = 1.f / std::max(m13*m23 - d * d, 1e-30f);
114 float a = x13.dot(x03), b = x23.dot(x03);
115 // the barycentric coordinates themselves
116 float w23 = invdet * (m23*a - d * b);
117 float w31 = invdet * (m13*b - d * a);
118 float w12 = 1 - w23 - w31;
119 if (w23 >= 0 && w31 >= 0 && w12 >= 0) { // if we're inside the triangle
120 return (x0 - (w23*x1 + w31 * x2 + w12 * x3)).norm();
121 }
122 else { // we have to clamp to one of the edges
123 if (w23 > 0) // this rules out edge 2-3 for us
124 return std::min(point_segment_distance(x0, x1, x2), point_segment_distance(x0, x1, x3));
125 else if (w31 > 0) // this rules out edge 1-3
126 return std::min(point_segment_distance(x0, x1, x2), point_segment_distance(x0, x2, x3));
127 else // w12 must be >0, ruling out edge 1-2
128 return std::min(point_segment_distance(x0, x1, x3), point_segment_distance(x0, x2, x3));
129 }
130 }
131
132 static void FSM_CheckNeighbour(const CArray<Vec3ui> &tri, const CArray<Vec3f> &vert,
133 CArray3f &phi, CArray3i &closest_tri,
134 const Vec3f &gx, int i0, int j0, int k0, int i1, int j1, int k1)
135 {
136 if (closest_tri(i1, j1, k1) >= 0) {
137 unsigned int p, q, r;
138 Vec3ui trijk = tri[closest_tri(i1, j1, k1)];
139 p = trijk[0];
140 q = trijk[1];
141 r = trijk[2];
142 float d = point_triangle_distance(gx, vert[p], vert[q], vert[r]);
143 if (d < phi(i0, j0, k0)) {
144 phi(i0, j0, k0) = d;
145 closest_tri(i0, j0, k0) = closest_tri(i1, j1, k1);
146 }
147 }
148 }
149
150 static void sweep(const CArray<Vec3ui> &tri, const CArray<Vec3f> &x,
151 CArray3f &phi, CArray3i &closest_tri, const Vec3f &origin, float dx,
152 int di, int dj, int dk)
153 {
154 int i0, i1;
155 if (di > 0) {
156 i0 = 1; i1 = phi.nx();
157 }
158 else {
159 i0 = phi.nx() - 2; i1 = -1;
160 }
161 int j0, j1;
162 if (dj > 0) {
163 j0 = 1; j1 = phi.ny();
164 }
165 else {
166 j0 = phi.ny() - 2; j1 = -1;
167 }
168 int k0, k1;
169 if (dk > 0) {
170 k0 = 1; k1 = phi.nz();
171 }
172 else {
173 k0 = phi.nz() - 2; k1 = -1;
174 }
175 for (int k = k0; k != k1; k += dk) for (int j = j0; j != j1; j += dj) for (int i = i0; i != i1; i += di) {
176 Vec3f gx(i*dx + origin[0], j*dx + origin[1], k*dx + origin[2]);
177 FSM_CheckNeighbour(tri, x, phi, closest_tri, gx, i, j, k, i - di, j, k);
178 FSM_CheckNeighbour(tri, x, phi, closest_tri, gx, i, j, k, i, j - dj, k);
179 FSM_CheckNeighbour(tri, x, phi, closest_tri, gx, i, j, k, i - di, j - dj, k);
180 FSM_CheckNeighbour(tri, x, phi, closest_tri, gx, i, j, k, i, j, k - dk);
181 FSM_CheckNeighbour(tri, x, phi, closest_tri, gx, i, j, k, i - di, j, k - dk);
182 FSM_CheckNeighbour(tri, x, phi, closest_tri, gx, i, j, k, i, j - dj, k - dk);
183 FSM_CheckNeighbour(tri, x, phi, closest_tri, gx, i, j, k, i - di, j - dj, k - dk);
184 }
185 }
186
187 // calculate twice signed area of triangle (0,0)-(x1,y1)-(x2,y2)
188 // return an SOS-determined sign (-1, +1, or 0 only if it's a truly degenerate triangle)
189 static int orientation(double x1, double y1, double x2, double y2, double &twice_signed_area)
190 {
191 twice_signed_area = y1 * x2 - x1 * y2;
192 if (twice_signed_area > 0) return 1;
193 else if (twice_signed_area < 0) return -1;
194 else if (y2 > y1) return 1;
195 else if (y2 < y1) return -1;
196 else if (x1 > x2) return 1;
197 else if (x1 < x2) return -1;
198 else return 0; // only true when x1==x2 and y1==y2
199 }
200
201 // robust test of (x0,y0) in the triangle (x1,y1)-(x2,y2)-(x3,y3)
202 // if true is returned, the barycentric coordinates are set in a,b,c.
203 static bool point_in_triangle_2d(double x0, double y0,
204 double x1, double y1, double x2, double y2, double x3, double y3,
205 double& a, double& b, double& c)
206 {
207 x1 -= x0; x2 -= x0; x3 -= x0;
208 y1 -= y0; y2 -= y0; y3 -= y0;
209 int signa = orientation(x2, y2, x3, y3, a);
210 if (signa == 0) return false;
211 int signb = orientation(x3, y3, x1, y1, b);
212 if (signb != signa) return false;
213 int signc = orientation(x1, y1, x2, y2, c);
214 if (signc != signa) return false;
215 double sum = a + b + c;
216 assert(sum != 0); // if the SOS signs match and are nonkero, there's no way all of a, b, and c are zero.
217 a /= sum;
218 b /= sum;
219 c /= sum;
220 return true;
221 }
222
223#define TRI_MIN(x, y, z) std::min(x, std::min(y, z))
224#define TRI_MAX(x, y, z) std::max(x, std::max(y, z))
225
226 template<typename TDataType>
228 {
229 const int exact_band = 1;
230 Real dx = this->varSpacing()->getValue();
231
232 phi.resize(ni, nj, nk);
233 phi.assign((ni + nj + nk)*dx);
234 //TODO: phi.assign((ni + nj + nk)*dx); // upper bound on distance
235 CArray3i closest_tri(ni, nj, nk);//CArray3i closest_tri(ni, nj, nk, -1);
236 closest_tri.assign(-1);
237 CArray3i intersection_count(ni, nj, nk); //CArray3i intersection_count(ni, nj, nk, 0); // intersection_count(i,j,k) is # of tri intersections in (i-1,i]x{j}x{k}
238 intersection_count.assign(0);
239 // we begin by initializing distances near the mesh, and figuring out intersection counts
240 for (uint t = 0; t < faceList.size(); ++t) {
241 uint p, q, r;
242 p = faceList[t][0];
243 q = faceList[t][1];
244 r = faceList[t][2];
245 // coordinates in grid to high precision
246 double fip = ((double)vertList[p][0] - origin[0]) / dx;
247 double fjp = ((double)vertList[p][1] - origin[1]) / dx;
248 double fkp = ((double)vertList[p][2] - origin[2]) / dx;
249 double fiq = ((double)vertList[q][0] - origin[0]) / dx;
250 double fjq = ((double)vertList[q][1] - origin[1]) / dx;
251 double fkq = ((double)vertList[q][2] - origin[2]) / dx;
252 double fir = ((double)vertList[r][0] - origin[0]) / dx;
253 double fjr = ((double)vertList[r][1] - origin[1]) / dx;
254 double fkr = ((double)vertList[r][2] - origin[2]) / dx;
255 // do distances nearby
256 int i0 = clamp(int(TRI_MIN(fip, fiq, fir)) - exact_band, 0, ni - 1);
257 int i1 = clamp(int(TRI_MAX(fip, fiq, fir)) + exact_band + 1, 0, ni - 1);
258 int j0 = clamp(int(TRI_MIN(fjp, fjq, fjr)) - exact_band, 0, nj - 1);
259 int j1 = clamp(int(TRI_MAX(fjp, fjq, fjr)) + exact_band + 1, 0, nj - 1);
260 int k0 = clamp(int(TRI_MIN(fkp, fkq, fkr)) - exact_band, 0, nk - 1);
261 int k1 = clamp(int(TRI_MAX(fkp, fkq, fkr)) + exact_band + 1, 0, nk - 1);
262 for (int k = k0; k <= k1; ++k) for (int j = j0; j <= j1; ++j) for (int i = i0; i <= i1; ++i) {
263 Vec3f gx(i*dx + origin[0], j*dx + origin[1], k*dx + origin[2]);
264 float d = point_triangle_distance(gx, vertList[p], vertList[q], vertList[r]);
265 if (d < phi(i, j, k)) {
266 phi(i, j, k) = d;
267 closest_tri(i, j, k) = t;
268 }
269 }
270 // and do intersection counts
271 j0 = clamp((int)std::ceil(TRI_MIN(fjp, fjq, fjr)), 0, nj - 1);
272 j1 = clamp((int)std::floor(TRI_MAX(fjp, fjq, fjr)), 0, nj - 1);
273 k0 = clamp((int)std::ceil(TRI_MIN(fkp, fkq, fkr)), 0, nk - 1);
274 k1 = clamp((int)std::floor(TRI_MAX(fkp, fkq, fkr)), 0, nk - 1);
275 for (int k = k0; k <= k1; ++k) for (int j = j0; j <= j1; ++j) {
276 double a, b, c;
277 if (point_in_triangle_2d(j, k, fjp, fkp, fjq, fkq, fjr, fkr, a, b, c)) {
278 double fi = a * fip + b * fiq + c * fir; // intersection i coordinate
279 int i_interval = int(std::ceil(fi)); // intersection is in (i_interval-1,i_interval]
280 if (i_interval < 0) ++intersection_count(0, j, k); // we enlarge the first interval to include everything to the -x direction
281 else if (i_interval < ni) ++intersection_count(i_interval, j, k);
282 // we ignore intersections that are beyond the +x side of the grid
283 }
284 }
285 }
286 // and now we fill in the rest of the distances with fast sweeping
287 for (unsigned int pass = 0; pass < 2; ++pass) {
288 sweep(faceList, vertList, phi, closest_tri, origin, dx, +1, +1, +1);
289 sweep(faceList, vertList, phi, closest_tri, origin, dx, -1, -1, -1);
290 sweep(faceList, vertList, phi, closest_tri, origin, dx, +1, +1, -1);
291 sweep(faceList, vertList, phi, closest_tri, origin, dx, -1, -1, +1);
292 sweep(faceList, vertList, phi, closest_tri, origin, dx, +1, -1, +1);
293 sweep(faceList, vertList, phi, closest_tri, origin, dx, -1, +1, -1);
294 sweep(faceList, vertList, phi, closest_tri, origin, dx, +1, -1, -1);
295 sweep(faceList, vertList, phi, closest_tri, origin, dx, -1, +1, +1);
296 }
297 // then figure out signs (inside/outside) from intersection counts
298 for (int k = 0; k < nk; ++k) for (int j = 0; j < nj; ++j) {
299 int total_count = 0;
300 for (int i = 0; i < ni; ++i) {
301 total_count += intersection_count(i, j, k);
302 if (total_count % 2 == 1) { // if parity of intersections so far is odd,
303 phi(i, j, k) = -phi(i, j, k); // we are inside the mesh
304 }
305 }
306 }
307 }
308
310}
#define TRI_MAX(x, y, z)
#define TRI_MIN(x, y, z)
#define DEFINE_CLASS(name)
Definition Object.h:140
#define IMPLEMENT_TCLASS(name, T1)
Definition Object.h:103
assert(queueCount >=1)
void setDistance(CArray3D< Real > distance)
void setSpace(const Coord p0, const Coord p1, Real h)
This is a CPU-based implementation of grid-based signed distance field (level set) generator for tria...
This is an implementation of AdditiveCCD based on peridyno.
Definition Array.h:25
static bool point_in_triangle_2d(double x0, double y0, double x1, double y1, double x2, double y2, double x3, double y3, double &a, double &b, double &c)
static float point_segment_distance(const Vec3f &x0, const Vec3f &x1, const Vec3f &x2)
DYN_FUNC T dot(Vector< T, 2 > const &U, Vector< T, 2 > const &V)
Definition SimpleMath.h:199
CArray3D< int > CArray3i
static int orientation(double x1, double y1, double x2, double y2, double &twice_signed_area)
static float point_triangle_distance(const Vec3f &x0, const Vec3f &x1, const Vec3f &x2, const Vec3f &x3)
static void FSM_CheckNeighbour(const CArray< Vec3ui > &tri, const CArray< Vec3f > &vert, CArray3f &phi, CArray3i &closest_tri, const Vec3f &gx, int i0, int j0, int k0, int i1, int j1, int k1)
static void sweep(const CArray< Vec3ui > &tri, const CArray< Vec3f > &x, CArray3f &phi, CArray3i &closest_tri, const Vec3f &origin, float dx, int di, int dj, int dk)
Vector< unsigned int, 3 > Vec3ui
Vector< float, 3 > Vec3f
Definition Vector3D.h:93
DYN_FUNC T clamp(const T &v, const T &lo, const T &hi)
Definition SimpleMath.h:42
Array< T, DeviceType::CPU > CArray
Definition Array.h:151
unsigned int uint
Definition VkProgram.h:14
CArray3D< float > CArray3f