// @(#)root/mathcore:$Id: Delaunay2D.h,v 1.00
// Author: Daniel Funke, Lorenzo Moneta

/*************************************************************************
 * Copyright (C) 2015 ROOT Math Team                                     *
 * All rights reserved.                                                  *
 *                                                                       *
 * For the licensing terms see $ROOTSYS/LICENSE.                         *
 * For the list of contributors see $ROOTSYS/README/CREDITS.             *
 *************************************************************************/

// Header file for class Delaunay2D

#ifndef ROOT_Math_Delaunay2D
#define ROOT_Math_Delaunay2D

//for testing purposes HAS_CGAL can be [un]defined here
//#define HAS_CGAL

//for testing purposes THREAD_SAFE can [un]defined here
//#define THREAD_SAFE


#ifndef ROOT_TNamed
#include "TNamed.h"
#endif

#include <map>
#include <vector>
#include <set>
#include <functional>

#ifdef HAS_CGAL
	/* CGAL uses the name PTR as member name in its Handle class
	 * but its a macro defined in mmalloc.h of ROOT
	 * Safe it, disable it and then re-enable it later on*/
	#pragma push_macro("PTR")
	#undef PTR

   #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
   #include <CGAL/Delaunay_triangulation_2.h>
   #include <CGAL/Triangulation_vertex_base_with_info_2.h>
   #include <CGAL/Interpolation_traits_2.h>
   #include <CGAL/natural_neighbor_coordinates_2.h>
   #include <CGAL/interpolation_functions.h>

	#pragma pop_macro("PTR")
#endif

#ifdef THREAD_SAFE
	#include<atomic> //atomic operations for thread safety
#endif


namespace ROOT {



   namespace Math {

/**

   Class to generate a Delaunay triangulation of a 2D set of points. 
   Algorithm based on **Triangle**, a two-dimensional quality mesh generator and 
   Delaunay triangulator from Jonathan Richard Shewchuk. 

   See [http://www.cs.cmu.edu/~quake/triangle.html]

   \ingroup MathCore
 */


class Delaunay2D  {

public:

	struct Triangle {
		double x[3];
		double y[3];
		UInt_t idx[3];

		#ifndef HAS_CGAL
		double invDenom; //see comment below in CGAL fall back section
		#endif
	};

	typedef std::vector<Triangle> Triangles;

public:

   
   Delaunay2D(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
   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);

   /// Return the Interpolated z value corresponding to the (x,y) point
   double  Interpolate(double x, double y);

   /// Find all triangles 
   void      FindAllTriangles();

   /// return the number of triangles 
   Int_t     NumberOfTriangles() const {return fNdt;}

   double  XMin() const {return fXNmin;}
   double  XMax() const {return fXNmax;}
   double  YMin() const {return fYNmin;}
   double  YMax() const {return fYNmax;}

   /// set z value to be returned for points  outside the region 
   void      SetZOuterValue(double z=0.) { fZout = z; }

   /// return the user defined Z-outer value
   double ZOuterValue() const { return fZout; } 

   // iterators on the found triangles 
   Triangles::const_iterator begin() const { return fTriangles.begin(); }
   Triangles::const_iterator end()  const { return fTriangles.end(); }


private:

   // internal methods

   
   inline double Linear_transform(double x, double offset, double factor){
	   return (x+offset)*factor;
   }

   /// internal function to normalize the points
   void DoNormalizePoints();

   /// internal function to find the triangle
   /// use Triangle or CGAL if flag is set 
   void DoFindTriangles();

   /// internal method to compute the interpolation
   double  DoInterpolateNormalized(double x, double y);


   
private:
   // class is not copyable
   Delaunay2D(const Delaunay2D&); // Not implemented
   Delaunay2D& operator=(const Delaunay2D&); // Not implemented

protected:

   //typedef std::function<double(double)>                  Transformer;

   Int_t       fNdt;         //!Number of Delaunay triangles found
   Int_t       fNpoints;     //!Number of data points 

   const double   *fX;           //!Pointer to X array (managed externally)
   const double   *fY;           //!Pointer to Y array 
   const double   *fZ;           //!Pointer to Z array 

   double    fXNmin;       //!Minimum value of fXN
   double    fXNmax;       //!Maximum value of fXN
   double    fYNmin;       //!Minimum value of fYN
   double    fYNmax;       //!Maximum value of fYN

   //Transformer xTransformer; //!transform x values to mapped space
   //Transformer yTransformer; //!transform y values to mapped space

   double    fOffsetX;      //!Normalization offset X
   double    fOffsetY;      //!Normalization offset Y

   double    fScaleFactorX; //!Normalization factor X
   double    fScaleFactorY; //!Normalization factor Y

   double    fZout;        //!Height for points lying outside the convex hull

#ifdef THREAD_SAFE

   enum class Initialization : char {UNINITIALIZED, INITIALIZING, INITIALIZED};
   std::atomic<Initialization> fInit; //!Indicate initialization state

#else
   Bool_t      fInit;        //!True if FindAllTriangels() has been performed
#endif


   Triangles   fTriangles;   //!Triangles of Triangulation

#ifdef HAS_CGAL

   //Functor class for accessing the function values/gradients
   	template< class PointWithInfoMap, typename ValueType >
   	struct Data_access : public std::unary_function< typename PointWithInfoMap::key_type,
   				 std::pair<ValueType, bool> >
   	{

   	  Data_access(const PointWithInfoMap& points, const ValueType * values)
   			  : _points(points), _values(values){};

   	  std::pair< ValueType, bool>
   	  operator()(const typename PointWithInfoMap::key_type& p) const {
   		typename PointWithInfoMap::const_iterator mit = _points.find(p);
   		if(mit!= _points.end())
   		  return std::make_pair(_values[mit->second], true);
   		return std::make_pair(ValueType(), false);
   	  };

   	  const PointWithInfoMap& _points;
   	  const ValueType * _values;
   	};

   	typedef CGAL::Exact_predicates_inexact_constructions_kernel  K;
   	typedef CGAL::Triangulation_vertex_base_with_info_2<uint, K> Vb;
   	typedef CGAL::Triangulation_data_structure_2<Vb>             Tds;
   	typedef CGAL::Delaunay_triangulation_2<K, Tds>               Delaunay;
   	typedef CGAL::Interpolation_traits_2<K>                      Traits;
   	typedef K::FT                                                Coord_type;
   	typedef K::Point_2                                           Point;
   	typedef std::map<Point, Vb::Info, K::Less_xy_2>              PointWithInfoMap;
   	typedef Data_access< PointWithInfoMap, double >            Value_access;

   Delaunay fCGALdelaunay; //! CGAL delaunay triangulation object
   PointWithInfoMap fNormalizedPoints; //! Normalized function values

#else // HAS_CGAL
   //fallback to triangle library

   /* Using barycentric coordinates for inTriangle test and interpolation
    *
    * Given triangle ABC and point P, P can be expressed by
    *
    * P.x = la * A.x + lb * B.x + lc * C.x
    * P.y = la * A.y + lb * B.y + lc * C.y
    *
    * with lc = 1 - la - lb
    *
    * P.x = la * A.x + lb * B.x + (1-la-lb) * C.x
    * P.y = la * A.y + lb * B.y + (1-la-lb) * C.y
    *
    * Rearranging yields
    *
    * la * (A.x - C.x) + lb * (B.x - C.x) = P.x - C.x
    * la * (A.y - C.y) + lb * (B.y - C.y) = P.y - C.y
    *
    * Thus
    *
    * la = ( (B.y - C.y)*(P.x - C.x) + (C.x - B.x)*(P.y - C.y) ) / ( (B.y - C.y)*(A.x - C.x) + (C.x - B.x)*(A.y - C.y) )
    * lb = ( (C.y - A.y)*(P.x - C.x) + (A.x - C.x)*(P.y - C.y) ) / ( (B.y - C.y)*(A.x - C.x) + (C.x - B.x)*(A.y - C.y) )
    * lc = 1 - la - lb
    *
    * We save the inverse denominator to speedup computation
    *
    * invDenom = 1 / ( (B.y - C.y)*(A.x - C.x) + (C.x - B.x)*(A.y - C.y) )
    *
    * P is in triangle (including edges if
    *
    * 0 <= [la, lb, lc] <= 1
    *
    * The interpolation of P.z is
    *
    * P.z = la * A.z + lb * B.z + lc * C.z
    *
    */

   std::vector<double> fXN; //! normalized X
   std::vector<double> fYN; //! normalized Y

   /* To speed up localisation of points a grid is layed over normalized space
    *
    * A reference to triangle ABC is added to _all_ grid cells that include ABC's bounding box
    */

   static const int fNCells = 25; //! number of cells to divide the normalized space
   double fXCellStep; //! inverse denominator to calculate X cell = fNCells / (fXNmax - fXNmin)
   double fYCellStep; //! inverse denominator to calculate X cell = fNCells / (fYNmax - fYNmin)
   std::set<UInt_t> fCells[(fNCells+1)*(fNCells+1)]; //! grid cells with containing triangles

   inline unsigned int Cell(UInt_t x, UInt_t y) const {
	   return x*(fNCells+1) + y;
   }

   inline int CellX(double x) const {
	   return (x - fXNmin) * fXCellStep;
   }

   inline int CellY(double y) const {
	   return (y - fYNmin) * fYCellStep;
   }

#endif //HAS_CGAL


};


} // namespace Math
} // namespace ROOT

      
#endif
 Delaunay2D.h:1
 Delaunay2D.h:2
 Delaunay2D.h:3
 Delaunay2D.h:4
 Delaunay2D.h:5
 Delaunay2D.h:6
 Delaunay2D.h:7
 Delaunay2D.h:8
 Delaunay2D.h:9
 Delaunay2D.h:10
 Delaunay2D.h:11
 Delaunay2D.h:12
 Delaunay2D.h:13
 Delaunay2D.h:14
 Delaunay2D.h:15
 Delaunay2D.h:16
 Delaunay2D.h:17
 Delaunay2D.h:18
 Delaunay2D.h:19
 Delaunay2D.h:20
 Delaunay2D.h:21
 Delaunay2D.h:22
 Delaunay2D.h:23
 Delaunay2D.h:24
 Delaunay2D.h:25
 Delaunay2D.h:26
 Delaunay2D.h:27
 Delaunay2D.h:28
 Delaunay2D.h:29
 Delaunay2D.h:30
 Delaunay2D.h:31
 Delaunay2D.h:32
 Delaunay2D.h:33
 Delaunay2D.h:34
 Delaunay2D.h:35
 Delaunay2D.h:36
 Delaunay2D.h:37
 Delaunay2D.h:38
 Delaunay2D.h:39
 Delaunay2D.h:40
 Delaunay2D.h:41
 Delaunay2D.h:42
 Delaunay2D.h:43
 Delaunay2D.h:44
 Delaunay2D.h:45
 Delaunay2D.h:46
 Delaunay2D.h:47
 Delaunay2D.h:48
 Delaunay2D.h:49
 Delaunay2D.h:50
 Delaunay2D.h:51
 Delaunay2D.h:52
 Delaunay2D.h:53
 Delaunay2D.h:54
 Delaunay2D.h:55
 Delaunay2D.h:56
 Delaunay2D.h:57
 Delaunay2D.h:58
 Delaunay2D.h:59
 Delaunay2D.h:60
 Delaunay2D.h:61
 Delaunay2D.h:62
 Delaunay2D.h:63
 Delaunay2D.h:64
 Delaunay2D.h:65
 Delaunay2D.h:66
 Delaunay2D.h:67
 Delaunay2D.h:68
 Delaunay2D.h:69
 Delaunay2D.h:70
 Delaunay2D.h:71
 Delaunay2D.h:72
 Delaunay2D.h:73
 Delaunay2D.h:74
 Delaunay2D.h:75
 Delaunay2D.h:76
 Delaunay2D.h:77
 Delaunay2D.h:78
 Delaunay2D.h:79
 Delaunay2D.h:80
 Delaunay2D.h:81
 Delaunay2D.h:82
 Delaunay2D.h:83
 Delaunay2D.h:84
 Delaunay2D.h:85
 Delaunay2D.h:86
 Delaunay2D.h:87
 Delaunay2D.h:88
 Delaunay2D.h:89
 Delaunay2D.h:90
 Delaunay2D.h:91
 Delaunay2D.h:92
 Delaunay2D.h:93
 Delaunay2D.h:94
 Delaunay2D.h:95
 Delaunay2D.h:96
 Delaunay2D.h:97
 Delaunay2D.h:98
 Delaunay2D.h:99
 Delaunay2D.h:100
 Delaunay2D.h:101
 Delaunay2D.h:102
 Delaunay2D.h:103
 Delaunay2D.h:104
 Delaunay2D.h:105
 Delaunay2D.h:106
 Delaunay2D.h:107
 Delaunay2D.h:108
 Delaunay2D.h:109
 Delaunay2D.h:110
 Delaunay2D.h:111
 Delaunay2D.h:112
 Delaunay2D.h:113
 Delaunay2D.h:114
 Delaunay2D.h:115
 Delaunay2D.h:116
 Delaunay2D.h:117
 Delaunay2D.h:118
 Delaunay2D.h:119
 Delaunay2D.h:120
 Delaunay2D.h:121
 Delaunay2D.h:122
 Delaunay2D.h:123
 Delaunay2D.h:124
 Delaunay2D.h:125
 Delaunay2D.h:126
 Delaunay2D.h:127
 Delaunay2D.h:128
 Delaunay2D.h:129
 Delaunay2D.h:130
 Delaunay2D.h:131
 Delaunay2D.h:132
 Delaunay2D.h:133
 Delaunay2D.h:134
 Delaunay2D.h:135
 Delaunay2D.h:136
 Delaunay2D.h:137
 Delaunay2D.h:138
 Delaunay2D.h:139
 Delaunay2D.h:140
 Delaunay2D.h:141
 Delaunay2D.h:142
 Delaunay2D.h:143
 Delaunay2D.h:144
 Delaunay2D.h:145
 Delaunay2D.h:146
 Delaunay2D.h:147
 Delaunay2D.h:148
 Delaunay2D.h:149
 Delaunay2D.h:150
 Delaunay2D.h:151
 Delaunay2D.h:152
 Delaunay2D.h:153
 Delaunay2D.h:154
 Delaunay2D.h:155
 Delaunay2D.h:156
 Delaunay2D.h:157
 Delaunay2D.h:158
 Delaunay2D.h:159
 Delaunay2D.h:160
 Delaunay2D.h:161
 Delaunay2D.h:162
 Delaunay2D.h:163
 Delaunay2D.h:164
 Delaunay2D.h:165
 Delaunay2D.h:166
 Delaunay2D.h:167
 Delaunay2D.h:168
 Delaunay2D.h:169
 Delaunay2D.h:170
 Delaunay2D.h:171
 Delaunay2D.h:172
 Delaunay2D.h:173
 Delaunay2D.h:174
 Delaunay2D.h:175
 Delaunay2D.h:176
 Delaunay2D.h:177
 Delaunay2D.h:178
 Delaunay2D.h:179
 Delaunay2D.h:180
 Delaunay2D.h:181
 Delaunay2D.h:182
 Delaunay2D.h:183
 Delaunay2D.h:184
 Delaunay2D.h:185
 Delaunay2D.h:186
 Delaunay2D.h:187
 Delaunay2D.h:188
 Delaunay2D.h:189
 Delaunay2D.h:190
 Delaunay2D.h:191
 Delaunay2D.h:192
 Delaunay2D.h:193
 Delaunay2D.h:194
 Delaunay2D.h:195
 Delaunay2D.h:196
 Delaunay2D.h:197
 Delaunay2D.h:198
 Delaunay2D.h:199
 Delaunay2D.h:200
 Delaunay2D.h:201
 Delaunay2D.h:202
 Delaunay2D.h:203
 Delaunay2D.h:204
 Delaunay2D.h:205
 Delaunay2D.h:206
 Delaunay2D.h:207
 Delaunay2D.h:208
 Delaunay2D.h:209
 Delaunay2D.h:210
 Delaunay2D.h:211
 Delaunay2D.h:212
 Delaunay2D.h:213
 Delaunay2D.h:214
 Delaunay2D.h:215
 Delaunay2D.h:216
 Delaunay2D.h:217
 Delaunay2D.h:218
 Delaunay2D.h:219
 Delaunay2D.h:220
 Delaunay2D.h:221
 Delaunay2D.h:222
 Delaunay2D.h:223
 Delaunay2D.h:224
 Delaunay2D.h:225
 Delaunay2D.h:226
 Delaunay2D.h:227
 Delaunay2D.h:228
 Delaunay2D.h:229
 Delaunay2D.h:230
 Delaunay2D.h:231
 Delaunay2D.h:232
 Delaunay2D.h:233
 Delaunay2D.h:234
 Delaunay2D.h:235
 Delaunay2D.h:236
 Delaunay2D.h:237
 Delaunay2D.h:238
 Delaunay2D.h:239
 Delaunay2D.h:240
 Delaunay2D.h:241
 Delaunay2D.h:242
 Delaunay2D.h:243
 Delaunay2D.h:244
 Delaunay2D.h:245
 Delaunay2D.h:246
 Delaunay2D.h:247
 Delaunay2D.h:248
 Delaunay2D.h:249
 Delaunay2D.h:250
 Delaunay2D.h:251
 Delaunay2D.h:252
 Delaunay2D.h:253
 Delaunay2D.h:254
 Delaunay2D.h:255
 Delaunay2D.h:256
 Delaunay2D.h:257
 Delaunay2D.h:258
 Delaunay2D.h:259
 Delaunay2D.h:260
 Delaunay2D.h:261
 Delaunay2D.h:262
 Delaunay2D.h:263
 Delaunay2D.h:264
 Delaunay2D.h:265
 Delaunay2D.h:266
 Delaunay2D.h:267
 Delaunay2D.h:268
 Delaunay2D.h:269
 Delaunay2D.h:270
 Delaunay2D.h:271
 Delaunay2D.h:272
 Delaunay2D.h:273
 Delaunay2D.h:274
 Delaunay2D.h:275
 Delaunay2D.h:276
 Delaunay2D.h:277
 Delaunay2D.h:278
 Delaunay2D.h:279
 Delaunay2D.h:280
 Delaunay2D.h:281
 Delaunay2D.h:282
 Delaunay2D.h:283
 Delaunay2D.h:284
 Delaunay2D.h:285
 Delaunay2D.h:286
 Delaunay2D.h:287
 Delaunay2D.h:288
 Delaunay2D.h:289
 Delaunay2D.h:290
 Delaunay2D.h:291
 Delaunay2D.h:292
 Delaunay2D.h:293
 Delaunay2D.h:294
 Delaunay2D.h:295
 Delaunay2D.h:296
 Delaunay2D.h:297
 Delaunay2D.h:298