64 fInit = Initialization::UNINITIALIZED;
69 if (
n == 0 || !
x || !
y || !z )
return;
72 xmin = *std::min_element(
x,
x+
n);
73 xmax = *std::max_element(
x,
x+
n);
75 ymin = *std::min_element(
y,
y+
n);
76 ymax = *std::max_element(
y,
y+
n);
133 if(
fInit.load(std::memory_order::memory_order_relaxed) == Initialization::INITIALIZED)
136 Initialization cState = Initialization::UNINITIALIZED;
137 if(
fInit.compare_exchange_strong(cState, Initialization::INITIALIZING,
138 std::memory_order::memory_order_release, std::memory_order::memory_order_relaxed))
160 fInit = Initialization::INITIALIZED;
161 }
else while(cState != Initialization::INITIALIZED) {
163 cState =
fInit.load(std::memory_order::memory_order_relaxed);
174void Delaunay2D::DonormalizePoints() {
180 fNormalizedPoints.insert(std::make_pair(p,
n));
186 fCGALdelaunay.insert(fNormalizedPoints.begin(), fNormalizedPoints.end());
188 std::transform(fCGALdelaunay.finite_faces_begin(),
189 fCGALdelaunay.finite_faces_end(), std::back_inserter(
fTriangles),
190 [] (
const Delaunay::Face face) -> Triangle {
194 auto transform = [&] (const unsigned int i) {
195 tri.x[i] = face.vertex(i)->point().x();
196 tri.y[i] = face.vertex(i)->point().y();
197 tri.idx[i] = face.vertex(i)->info();
210double Delaunay2D::DoInterpolateNormalized(
double xx,
double yy)
222 std::vector<std::pair<Point, Coord_type> > coords;
223 auto nn = CGAL::natural_neighbor_coordinates_2(fCGALdelaunay, p,
224 std::back_inserter(coords));
233 Coord_type res = CGAL::linear_interpolation(coords.begin(), coords.end(),
234 nn.second, Value_access(fNormalizedPoints, fZ));
244void Delaunay2D::DoNormalizePoints() {
245 for (
Int_t n = 0;
n < fNpoints;
n++) {
246 fXN.push_back(Linear_transform(fX[
n], fOffsetX, fScaleFactorX));
247 fYN.push_back(Linear_transform(fY[
n], fOffsetY, fScaleFactorY));
251 fXCellStep = fNCells / (fXNmax - fXNmin);
252 fYCellStep = fNCells / (fYNmax - fYNmin);
256void Delaunay2D::DoFindTriangles() {
259 s.pointlist =
nullptr;
260 s.pointattributelist =
nullptr;
261 s.pointmarkerlist =
nullptr;
262 s.numberofpoints = 0;
263 s.numberofpointattributes = 0;
265 s.trianglelist =
nullptr;
266 s.triangleattributelist =
nullptr;
267 s.trianglearealist =
nullptr;
268 s.neighborlist =
nullptr;
269 s.numberoftriangles = 0;
270 s.numberofcorners = 0;
271 s.numberoftriangleattributes = 0;
273 s.segmentlist =
nullptr;
274 s.segmentmarkerlist =
nullptr;
275 s.numberofsegments = 0;
277 s.holelist =
nullptr;
280 s.regionlist =
nullptr;
281 s.numberofregions = 0;
283 s.edgelist =
nullptr;
284 s.edgemarkerlist =
nullptr;
285 s.normlist =
nullptr;
290 if(s.pointlist !=
nullptr)
free(s.pointlist);
291 if(s.pointattributelist !=
nullptr)
free(s.pointattributelist);
292 if(s.pointmarkerlist !=
nullptr)
free(s.pointmarkerlist);
294 if(s.trianglelist !=
nullptr)
free(s.trianglelist);
295 if(s.triangleattributelist !=
nullptr)
free(s.triangleattributelist);
296 if(s.trianglearealist !=
nullptr)
free(s.trianglearealist);
297 if(s.neighborlist !=
nullptr)
free(s.neighborlist);
299 if(s.segmentlist !=
nullptr)
free(s.segmentlist);
300 if(s.segmentmarkerlist !=
nullptr)
free(s.segmentmarkerlist);
302 if(s.holelist !=
nullptr)
free(s.holelist);
304 if(s.regionlist !=
nullptr)
free(s.regionlist);
306 if(s.edgelist !=
nullptr)
free(s.edgelist);
307 if(s.edgemarkerlist !=
nullptr)
free(s.edgemarkerlist);
308 if(s.normlist !=
nullptr)
free(s.normlist);
312 initStruct(in); initStruct(out);
319 for (
Int_t i = 0; i < fNpoints; ++i) {
330 auto transform = [&] (
const unsigned int v) {
351 tri.
invDenom = 1 / ( (tri.
y[1] - tri.
y[2])*(tri.
x[0] - tri.
x[2]) + (tri.
x[2] - tri.
x[1])*(tri.
y[0] - tri.
y[2]) );
355 auto bx = std::minmax({tri.
x[0], tri.
x[1], tri.
x[2]});
356 auto by = std::minmax({tri.
y[0], tri.
y[1], tri.
y[2]});
358 unsigned int cellXmin = CellX(bx.first);
359 unsigned int cellXmax = CellX(bx.second);
361 unsigned int cellYmin = CellY(by.first);
362 unsigned int cellYmax = CellY(by.second);
364 for(
unsigned int i = cellXmin; i <= cellXmax; ++i) {
365 for(
unsigned int j = cellYmin; j <= cellYmax; ++j) {
367 fCells[Cell(i,j)].insert(t);
372 freeStruct(in); freeStruct(out);
379double Delaunay2D::DoInterpolateNormalized(
double xx,
double yy)
386 auto bayCoords = [&] (
const unsigned int t) -> std::tuple<double, double, double> {
387 double la = ( (fTriangles[t].y[1] - fTriangles[t].y[2])*(xx - fTriangles[t].
x[2])
388 + (fTriangles[t].x[2] - fTriangles[t].x[1])*(yy - fTriangles[t].
y[2]) ) * fTriangles[t].invDenom;
389 double lb = ( (fTriangles[t].y[2] - fTriangles[t].y[0])*(xx - fTriangles[t].
x[2])
390 + (fTriangles[t].x[0] - fTriangles[t].x[2])*(yy - fTriangles[t].
y[2]) ) * fTriangles[t].invDenom;
392 return std::make_tuple(la, lb, (1 - la - lb));
395 auto inTriangle = [] (
const std::tuple<double, double, double> & coords) ->
bool {
396 return std::get<0>(coords) >= 0 && std::get<1>(coords) >= 0 && std::get<2>(coords) >= 0;
402 if(cX < 0 || cX > fNCells || cY < 0 || cY > fNCells)
405 for(
unsigned int t : fCells[Cell(cX, cY)]){
406 auto coords = bayCoords(t);
408 if(inTriangle(coords)){
410 return std::get<0>(coords) * fZ[fTriangles[t].idx[0]]
411 + std::get<1>(coords) * fZ[fTriangles[t].idx[1]]
412 + std::get<2>(coords) * fZ[fTriangles[t].idx[2]];
double fXNmin
Pointer to Z array.
double DoInterpolateNormalized(double x, double y)
internal method to compute the interpolation
double Interpolate(double x, double y)
Return the Interpolated z value corresponding to the (x,y) point.
const double * fZ
Pointer to Y array.
double fXNmax
Minimum value of fXN.
double fZout
Normalization factor Y.
double fOffsetY
Normalization offset X.
Triangles fTriangles
True if FindAllTriangels() has been performed.
double Linear_transform(double x, double offset, double factor)
double fYCellStep
inverse denominator to calculate X cell = fNCells / (fXNmax - fXNmin)
Delaunay2D(int n, const double *x, const double *y, const double *z, double xmin=0, double xmax=0, double ymin=0, double ymax=0)
class constructor from array of data points
void FindAllTriangles()
Find all triangles.
double fXCellStep
number of cells to divide the normalized space
double fYNmin
Maximum value of fXN.
double fOffsetX
Maximum value of fYN.
Int_t fNpoints
Number of Delaunay triangles found.
Bool_t fInit
Height for points lying outside the convex hull.
const double * fY
Pointer to X array (managed externally)
void SetInputPoints(int n, const double *x, const double *y, const double *z, double xmin=0, double xmax=0, double ymin=0, double ymax=0)
set the input points for building the graph
const double * fX
Number of data points.
void DoNormalizePoints()
internal function to normalize the points
double fYNmax
Minimum value of fYN.
double fScaleFactorY
Normalization factor X.
double fScaleFactorX
Normalization offset Y.
void DoFindTriangles()
internal function to find the triangle use Triangle or CGAL if flag is set
Namespace for new Math classes and functions.
tbb::task_arena is an alias of tbb::interface7::task_arena, which doesn't allow to forward declare tb...
void triangulate(char *triswitches, struct triangulateio *in, struct triangulateio *out, struct triangulateio *vorout)