libtin/tin.h
2021-09-21 22:26:58 +08:00

463 lines
14 KiB
C++
Raw Blame History

This file contains ambiguous Unicode characters

This file contains Unicode characters that might be confused with other characters. If you think that this is intentional, you can safely ignore this warning. Use the Escape button to reveal them.

/**
* @defgroup TIN
*
* @brief Generation of a Triangular Irregular Network (TIN) from a dense DEM grid.
*
* @author Yi Zhang
* @date 2021-09-16
*/
#ifndef _TIN_DELAUNAY_H
#define _TIN_DELAUNAY_H
#include "cmath"
#include "vector"
#include "algorithm"
#include "iostream"
#define ZERO 1e-5
// Start vertex definition
struct vertex2dc
{
unsigned int id; // index of the vertex
double x, y; // position of the vertex
double elev; // elevation at the vertex
vertex2dc() : x(NAN), y(NAN), elev(NAN), id(0) {}
vertex2dc(double inx, double iny, double inelev, unsigned int inid) {set(inx, iny, inelev, inid);}
void set(double inx, double iny, double inelev, unsigned int inid)
{
x = inx; y = iny; elev = inelev; id = inid;
return;
}
};
bool operator ==(const vertex2dc &a, const vertex2dc &b) // overload the == operator for vertex2dc type
{
if(fabs(a.x - b.x) <= ZERO && fabs(a.y - b.y) <= ZERO)
{
return true;
}
return false;
}
bool is_collinear(vertex2dc *a_ptr, vertex2dc *b_ptr, vertex2dc *c_ptr) // Test if the three points are on the same line
{
// (y3y1)(x2x1)(y2y1)(x3x1)
if (fabs((c_ptr->y - a_ptr->y)*(b_ptr->x - a_ptr->x) - (b_ptr->y - a_ptr->y)*(c_ptr->x - a_ptr->x)) <= ZERO)
{
return true;
}
return false;
}
// End vertex definition
// Start edge definition
struct edge
{
vertex2dc *vert[2]; // vertex of the edge
edge() {vert[0] = vert[1] = nullptr;}
edge(vertex2dc *v0ptr, vertex2dc *v1ptr) {set(v0ptr, v1ptr);}
void set(vertex2dc *v0ptr, vertex2dc *v1ptr)
{
vert[0] = v0ptr; vert[1] = v1ptr;
return;
}
};
bool operator ==(const edge &a, const edge &b) // overload the == operator for edge type
{
if((a.vert[0] == b.vert[0] && a.vert[1] == b.vert[1]) ||
(a.vert[0] == b.vert[1] && a.vert[1] == b.vert[0]))
{
return true;
}
return false;
}
// End edge definition
// Start triangle definition
struct dem_point;
struct triangle
{
vertex2dc *vert[3]; // vertex of the triangle
double cx, cy; // center of the triangle's circumcircle
double cr; // radius of the circumcircle
std::vector<dem_point*> circum_dem;
triangle() {vert[0] = vert[1] = vert[2] = nullptr;}
triangle(vertex2dc *v0ptr, vertex2dc *v1ptr, vertex2dc *v2ptr) {set(v0ptr, v1ptr, v2ptr);}
void set(vertex2dc *v0ptr, vertex2dc *v1ptr, vertex2dc *v2ptr)
{
vert[0] = v0ptr; vert[1] = v1ptr; vert[2] = v2ptr;
double s = 0.5 / ((vert[1]->x - vert[0]->x) * (vert[2]->y - vert[0]->y) - (vert[1]->y - vert[0]->y) * (vert[2]->x - vert[0]->x));
double m = vert[1]->x * vert[1]->x - vert[0]->x * vert[0]->x + vert[1]->y * vert[1]->y - vert[0]->y * vert[0]->y;
double u = vert[2]->x * vert[2]->x - vert[0]->x * vert[0]->x + vert[2]->y * vert[2]->y - vert[0]->y * vert[0]->y;
cx = ((vert[2]->y - vert[0]->y) * m + (vert[0]->y - vert[1]->y) * u) * s;
cy = ((vert[0]->x - vert[2]->x) * m + (vert[1]->x - vert[0]->x) * u) * s;
cr = (vert[0]->x - cx) * (vert[0]->x - cx) + (vert[0]->y - cy) * (vert[0]->y - cy); // not need to sqrt() here
return;
}
bool bound_location(double inx, double iny) // Test if the location is inside the triangle
{
double l1x, l1y, l2x, l2y;
for (int i = 0; i < 3; i++)
{
l1x = vert[(i+1)%3]->x - vert[i]->x;
l1y = vert[(i+1)%3]->y - vert[i]->y;
l2x = inx - vert[i]->x;
l2y = iny - vert[i]->y;
if ((l1x*l2y - l1y*l2x) < 0) // This condition includes points on the triangle's edge
{
return false;
}
}
return true;
}
double interpolate(double inx, double iny) // Interpolate the elevation of the given location inside the triangle
{
double a1 = 0.5 * ((vert[1]->x - inx)*(vert[2]->y - iny) - (vert[1]->y - iny)*(vert[2]->x - inx));
double a2 = 0.5 * ((vert[2]->x - inx)*(vert[0]->y - iny) - (vert[2]->y - iny)*(vert[0]->x - inx));
double a3 = 0.5 * ((vert[0]->x - inx)*(vert[1]->y - iny) - (vert[0]->y - iny)*(vert[1]->x - inx));
return (a1*vert[0]->elev + a2*vert[1]->elev + a3*vert[2]->elev)/(a1 + a2 + a3);
}
};
// End triangle definition
// Start DEM definition
struct dem_point
{
double x, y; // position of the DEM location
double elev; // elevation at the DEM location
double err; // error of the TIN with respect to the elevation
triangle *host; // host triangle of the DEM location
std::vector<triangle*> circum_host; // triangles which circumcircles include the location
dem_point() : x(NAN), y(NAN), elev(NAN), host(nullptr) {}
dem_point(double inx, double iny, double inelev) {set(inx, iny, inelev);}
void set(double inx, double iny, double inelev)
{
x = inx; y = iny; elev = inelev; host = nullptr;
return;
}
};
bool compare_dem_point(dem_point *a, dem_point *b)
{
if (a->err > b->err) return true;
return false;
}
// End DEM definition
/**
* @brief Generate the TIN from the DEM grid
*
* @param[in] dem Input DEM grid (Ordered from lower left corner to the upper right corner)
* @param[in] xmin The minimal coordinate of the DEM grid on the x-axis
* @param[in] xmax The maximal coordinate of the DEM grid on the x-axis
* @param[in] ymin The minimal coordinate of the DEM grid on the y-axis
* @param[in] ymax The maximal coordinate of the DEM grid on the y-axis
* @param[in] dx Data spacing of the DEM grid on the x-axis
* @param[in] dy Data spacing of the DEM grid on the y-axis
* @param out_verts The output vector of vertex's pointers. The user need to destroy the memories allocated by the function before destroy the vector
* @param out_tris The output vector of triangle's pointers. The user need to destroy the memories allocated by the function before destroy the vector
* @param[in] maxi_err Threshold to quit the algorithm. The default is 1e-0
* @param[in] err_records If this pointer is not NULL, record maximal error values after each insertion of vertex.
*/
void dem2tin(const std::vector<double> &dem, double xmin, double xmax, double ymin, double ymax,
double dx, double dy, std::vector<vertex2dc*> &out_verts, std::vector<triangle*> &out_tris,
double maxi_err = 1e-0, std::vector<double> *err_records = nullptr)
{
if (!out_verts.empty()) out_verts.clear();
if (!out_tris.empty()) out_tris.clear();
if (err_records != nullptr && !err_records->empty()) err_records->clear();
if (dx <= 0.0 || dy <= 0.0 || maxi_err <= 0.0) return;
if (xmin >= xmax || ymin >= ymax || (xmin + dx) > xmax || (ymin + dy) > ymax) return;
int xnum = round((xmax - xmin)/dx) + 1;
int ynum = round((ymax - ymin)/dy) + 1;
if (dem.size() != xnum*ynum) return;
// Prepare the DEM points
dem_point *tmp_dem = nullptr;;
std::vector<dem_point*> dem_grid(xnum*ynum);
std::vector<dem_point*>::iterator d_iter;
for (int i = 0; i < ynum; ++i)
{
for (int j = 0; j < xnum; ++j)
{
dem_grid[j + i*xnum] = new dem_point(xmin + dx*j, ymin + dy*i, dem[j + i*xnum]);
}
}
vertex2dc *tmp_vert = nullptr;
tmp_vert = new vertex2dc(xmin, ymin, dem_grid[0]->elev, out_verts.size()); // lower left corner
out_verts.push_back(tmp_vert);
d_iter = dem_grid.begin();
tmp_dem = *d_iter; delete tmp_dem;
dem_grid.erase(d_iter);
tmp_vert = new vertex2dc(xmax, ymin, dem_grid[xnum-2]->elev, out_verts.size()); // lower right corner. Note the first location is already erased
out_verts.push_back(tmp_vert);
d_iter = dem_grid.begin() + (xnum - 2);
tmp_dem = *d_iter; delete tmp_dem;
dem_grid.erase(d_iter);
tmp_vert = new vertex2dc(xmax, ymax, dem_grid[xnum*ynum-3]->elev, out_verts.size()); // upper right corner. Note the first two locations are already erased
out_verts.push_back(tmp_vert);
d_iter = dem_grid.begin() + (xnum*ynum - 3);
tmp_dem = *d_iter; delete tmp_dem;
dem_grid.erase(d_iter);
tmp_vert = new vertex2dc(xmin, ymax, dem_grid[xnum*(ynum-1) - 2]->elev, out_verts.size()); // upper left corner. Note the first two locations are already erased
out_verts.push_back(tmp_vert);
d_iter = dem_grid.begin() + (xnum*(ynum-1) - 2);
tmp_dem = *d_iter; delete tmp_dem;
dem_grid.erase(d_iter);
triangle *tmp_tri = nullptr;
std::vector<triangle*> cnst_tri, new_tri;
std::vector<triangle*>::iterator t_iter;
if (!is_collinear(out_verts[0], out_verts[1], out_verts[2])) // Do not create triangle if the vertexes are collinear
{
tmp_tri = new triangle(out_verts[0], out_verts[1], out_verts[2]); // order the vertex anti-clock wise
out_tris.push_back(tmp_tri);
}
if (!is_collinear(out_verts[0], out_verts[2], out_verts[3]))
{
tmp_tri = new triangle(out_verts[0], out_verts[2], out_verts[3]); // order the vertex anti-clock wise
out_tris.push_back(tmp_tri);
}
// Find host triangle for all DEM locations
for (int i = 0; i < dem_grid.size(); ++i)
{
for (int t = 0; t < out_tris.size(); ++t)
{
if (out_tris[t]->bound_location(dem_grid[i]->x, dem_grid[i]->y))
{
dem_grid[i]->host = out_tris[t];
break; // already found, no need to search more
}
}
}
// Find circum_host triangles for all DEM locations
double dist;
for (int i = 0; i < dem_grid.size(); ++i)
{
for (int t = 0; t < out_tris.size(); ++t)
{
dist = (out_tris[t]->cx - dem_grid[i]->x) * (out_tris[t]->cx - dem_grid[i]->x)
+ (out_tris[t]->cy - dem_grid[i]->y) * (out_tris[t]->cy - dem_grid[i]->y);
if ((dist - out_tris[t]->cr) <= ZERO) // Points on the circumcircle are also included
{
dem_grid[i]->circum_host.push_back(out_tris[t]);
out_tris[t]->circum_dem.push_back(dem_grid[i]);
// no beak here. There might be more than one triangle's circumcircle includes the DEM location
}
}
}
// loop all DEM data to find the location with maximal error
for (int i = 0; i < dem_grid.size(); ++i)
{
dem_grid[i]->err = fabs(dem_grid[i]->host->interpolate(dem_grid[i]->x, dem_grid[i]->y) - dem_grid[i]->elev);
}
// Sort dem_grid in the desceding order with respect to the error
std::sort(dem_grid.begin(), dem_grid.end(), compare_dem_point);
bool removed;
edge tmp_edge;
std::vector<edge> cnst_edge;
std::vector<edge>::iterator e_iter;
while (dem_grid[0]->err >= maxi_err) // quit til the threshold is meet
{
if (err_records != nullptr)
{
err_records->push_back(dem_grid[0]->err);
}
// create a new vertex
tmp_vert = new vertex2dc(dem_grid[0]->x, dem_grid[0]->y, dem_grid[0]->elev, out_verts.size());
out_verts.push_back(tmp_vert);
// Move triangles which circumcircles include the new vertex to the cnst_tri and remove it from out_tris
cnst_tri.clear();
for (int i = 0; i < dem_grid[0]->circum_host.size(); ++i)
{
cnst_tri.push_back(dem_grid[0]->circum_host[i]);
}
for (int c = 0; c < cnst_tri.size(); ++c)
{
for (t_iter = out_tris.begin(); t_iter != out_tris.end(); )
{
tmp_tri = *t_iter;
if (cnst_tri[c] == tmp_tri)
{
t_iter = out_tris.erase(t_iter);
break; // no need to search more
}
else t_iter++;
}
}
// remove cnst_tri from its circumed DEM's circum triangle list
for (int c = 0; c < cnst_tri.size(); ++c)
{
for (int i = 0; i < cnst_tri[c]->circum_dem.size(); ++i)
{
tmp_dem = cnst_tri[c]->circum_dem[i];
for (t_iter = tmp_dem->circum_host.begin(); t_iter != tmp_dem->circum_host.end(); )
{
if (cnst_tri[c] == *t_iter)
{
t_iter = tmp_dem->circum_host.erase(t_iter);
break;
}
else t_iter++;
}
}
}
// remove dem_grid[0] from its circumed triangle's circum DEM list
for (int c = 0; c < cnst_tri.size(); ++c)
{
for (d_iter = cnst_tri[c]->circum_dem.begin(); d_iter != cnst_tri[c]->circum_dem.end(); )
{
if (dem_grid[0] == *d_iter)
{
d_iter = cnst_tri[c]->circum_dem.erase(d_iter);
break;
}
else d_iter++;
}
}
// clear host and circumcircle triangles for the used DEM location
d_iter = dem_grid.begin();
tmp_dem = *d_iter; tmp_dem->circum_host.clear(); delete tmp_dem;
dem_grid.erase(d_iter);
// loop to remove duplicate edges
cnst_edge.clear();
for (int c = 0; c < cnst_tri.size(); ++c)
{
for (int e = 0; e < 3; ++e)
{
tmp_edge.set(cnst_tri[c]->vert[e], cnst_tri[c]->vert[(e+1)%3]);
removed = false;
for (e_iter = cnst_edge.begin(); e_iter != cnst_edge.end(); )
{
if (tmp_edge == *e_iter) // duplicate edge, remove from cnst_edge
{
e_iter = cnst_edge.erase(e_iter);
removed = true;
break; // no need to search more
}
else e_iter++;
}
if (!removed) // not a duplicate edge, add to the cnst_edge
{
cnst_edge.push_back(tmp_edge);
}
}
}
// construct new triangles and add to out_tris
new_tri.clear();
for (int c = 0; c < cnst_edge.size(); ++c)
{
if (!is_collinear(cnst_edge[c].vert[0], cnst_edge[c].vert[1], tmp_vert)) // Do not create triangle if the vertexes are collinear
{
tmp_tri = new triangle(cnst_edge[c].vert[0], cnst_edge[c].vert[1], tmp_vert); // order the vertex anti-clock wise
out_tris.push_back(tmp_tri);
new_tri.push_back(tmp_tri);
}
}
// loop all DEM data to update host triangles
for (int c = 0; c < cnst_tri.size(); ++c)
{
for (int i = 0; i < cnst_tri[c]->circum_dem.size(); ++i)
{
tmp_dem = cnst_tri[c]->circum_dem[i];
for (int n = 0; n < new_tri.size(); ++n) // search in newly created triangles to find new host
{
if (new_tri[n]->bound_location(tmp_dem->x, tmp_dem->y))
{
tmp_dem->host = new_tri[n];
tmp_dem->err = fabs(new_tri[n]->interpolate(tmp_dem->x, tmp_dem->y) - tmp_dem->elev);
break; // already found, no need to search more
}
}
}
}
// Find circum_host triangles for all DEM locations
// cnst_tri's circum area doesn't over cover new_tri's circum area
for (int i = 0; i < dem_grid.size(); ++i)
{
for (int n = 0; n < new_tri.size(); ++n) // search in newly created triangles to find new circumcircle triangles
{
dist = (new_tri[n]->cx - dem_grid[i]->x) * (new_tri[n]->cx - dem_grid[i]->x)
+ (new_tri[n]->cy - dem_grid[i]->y) * (new_tri[n]->cy - dem_grid[i]->y);
if ((dist - new_tri[n]->cr) <= ZERO) // Points on the circumcircle are also included
{
new_tri[n]->circum_dem.push_back(dem_grid[i]);
dem_grid[i]->circum_host.push_back(new_tri[n]);
// no beak here. There might be more than one triangle's circumcircle includes the DEM location
}
}
}
// destroy memories used by cnst_edge
for (int c = 0; c < cnst_tri.size(); ++c)
{
tmp_tri = cnst_tri[c];
tmp_tri->circum_dem.clear();
delete tmp_tri; tmp_tri = nullptr;
}
// Sort dem_grid in the desceding order with respect to the error
std::sort(dem_grid.begin(), dem_grid.end(), compare_dem_point);
}
if (err_records != nullptr)
{
err_records->push_back(dem_grid[0]->err);
}
// destroy remaining DEM data
for (int i = 0; i < dem_grid.size(); ++i)
{
tmp_dem = dem_grid[i];
delete tmp_dem; tmp_dem = nullptr;
}
return;
}
#endif // _TIN_DELAUNAY_H