Geogram Version 1.9.6-rc
A programming library of geometric algorithms
Loading...
Searching...
No Matches
GEO::OptimalTransportMap2d Class Reference

Computes semi-discrete optimal transport maps. More...

#include <exploragram/optimal_transport/optimal_transport_2d.h>

Inheritance diagram for GEO::OptimalTransportMap2d:
GEO::OptimalTransportMap

Public Member Functions

 OptimalTransportMap2d (Mesh *mesh, const std::string &delaunay="BPOW2d", bool BRIO=false)
 OptimalTransportMap2d constructor.
 
 ~OptimalTransportMap2d () override
 OptimalTransportMap destructor.
 
void get_RVD (Mesh &M) override
 Computes a mesh with the restricted Voronoi diagram.
 
void compute_Laguerre_centroids (double *centroids) override
 Computes the centroids of the Laguerre cells.
 
double total_mesh_mass () const
 Gets the total mass of the mesh.
 
- Public Member Functions inherited from GEO::OptimalTransportMap
 OptimalTransportMap (index_t dimension, Mesh *mesh, const std::string &delaunay="default", bool BRIO=false)
 OptimalTransportMap constructor.
 
virtual ~OptimalTransportMap ()
 OptimalTransportMap destructor.
 
index_t dimension () const
 Gets the dimension.
 
Meshmesh ()
 Gets the mesh.
 
void set_Newton (bool x)
 Sets whether Newton algorithm should be used.
 
void set_points (index_t nb_points, const double *points, index_t stride=0)
 Sets the points that define the target distribution.
 
void set_air_particles (index_t nb_air_particles, const double *air_particles, index_t stride, double air_fraction)
 Sets the air particles that define the volume occupied by the free space.
 
double air_fraction () const
 Gets the air fraction.
 
index_t nb_air_particles () const
 Gets the number of air particles.
 
void set_nu (index_t i, double nu)
 Sets the desired mass at one of the Diracs.
 
void set_Laguerre_centroids (double *x)
 Specifies a user vector where the centroids of the Laguerre cells will be stored after computing transport.
 
void set_epsilon (double eps)
 Sets the maximum error.
 
void set_linsolve_epsilon (double eps)
 Sets the tolerance for linear solve.
 
void set_linsolve_maxiter (index_t maxiter)
 Sets the maximum number of iterations for linear solve.
 
void set_linesearch_maxiter (index_t maxiter)
 Sets the maximum number of line search iterations.
 
void set_linesearch_init_iter (index_t init_iter)
 Sets the number of steplength reductions to be done at the first iteration.
 
void set_regularization (double eps_reg)
 Sets the weight of the regularization term.
 
void set_linear_solver (OTLinearSolver solver)
 Specifies whether a direct solver should be used.
 
void optimize (index_t max_iterations)
 Computes the weights that realize the optimal transport map between the source mesh and the target pointset.
 
void set_verbose (bool x)
 Enable/disable messages during optimization.
 
void optimize_full_Newton (index_t max_iterations, index_t n=0)
 Computes the weights that realize the optimal transport map between the source mesh and the target pointset.
 
void optimize_level (index_t b, index_t e, index_t max_iterations)
 Optimizes one level of the multilevel algorithm.
 
void optimize_levels (const vector< index_t > &levels, index_t max_iterations)
 Multi-level optimization.
 
index_t nb_points () const
 Gets the number of points.
 
const double * point_ptr (index_t i) const
 Gets a point.
 
double weight (index_t i) const
 Gets weight of a point.
 
void set_weight (index_t i, double val)
 Sets a weight of a point.
 
double potential (index_t i) const
 Gets the d+1-th coordinate of the embedding for a point.
 
RestrictedVoronoiDiagramRVD ()
 Gets the restricted Voronoi diagram.
 
void set_save_RVD_iter (bool x, bool show_RVD_seed=false, bool last_iter_only=false)
 Sets whether the restricted Voronoi diagram at each iteration should be saved.
 
void update_sparsity_pattern ()
 Updates the sparsity pattern of the Hessian right after a new Laguerre diagram was computed.
 
void new_linear_system (index_t n, double *x)
 Starts a new linear system.
 
void add_ij_coefficient (index_t i, index_t j, double a)
 Adds a coefficient to the matrix of the system.
 
void add_i_right_hand_side (index_t i, double a)
 Adds a coefficient to the right hand side.
 
void solve_linear_system ()
 Solves a linear system.
 
void set_initial_weight (index_t i, double w)
 Sets the initial value of the weight associated with one of the points.
 
double total_mass () const
 Gets the total mass of the domain.
 
void compute_P1_Laplacian (const double *weights, NLMatrix Laplacian, double *measures)
 Computes the P1 Laplacian of the Laguerre cells.
 

Public Attributes

GEOGen::Polygon work_
 Used by clipping operations.
 
GEOGen::Polygon clipped_
 Used by clipping operations.
 

Protected Member Functions

void call_callback_on_RVD () override
 Calls the callback for each intersection between a Laguerre cell and a simplex of the background mesh.
 
- Protected Member Functions inherited from GEO::OptimalTransportMap
double nu (index_t p) const
 Gets the mass of the Dirac associated with point p.
 
virtual void newiteration ()
 Callback for the numerical solver.
 
void save_RVD (index_t id)
 Saves the RVD at each iteration if specified on command line (just for debugging/ explaining the algorithm).
 
void funcgrad (index_t n, double *w, double &f, double *g)
 Computes the objective function and its gradient.
 
void eval_func_grad_Hessian (index_t n, const double *w, double &f, double *g)
 Computes the objective function, its gradient and its Hessian.
 
double gradient_threshold (index_t n) const
 Computes the stopping criterion of the solver.
 

Additional Inherited Members

- Static Public Member Functions inherited from GEO::OptimalTransportMap
static void funcgrad_CB (index_t n, double *x, double &f, double *g)
 Callback for the numerical solver.
 
static void newiteration_CB (index_t n, const double *x, double f, const double *g, double gnorm)
 Callback for the numerical solver.
 
- Protected Attributes inherited from GEO::OptimalTransportMap
index_t dimension_
 
index_t dimp1_
 dimension_ + 1
 
Meshmesh_
 
Delaunay_var delaunay_
 
RestrictedVoronoiDiagram_var RVD_
 
vector< double > points_dimp1_
 
vector< double > weights_
 
double total_mass_
 
double constant_nu_
 Value of one of the Diracs if cte.
 
vector< double > nu_
 Value of all the Diracs.
 
double epsilon_
 Acceptable relative deviation for the measure of a cell.
 
index_t current_call_iter_
 
Callbackcallback_
 
std::string last_stats_
 
bool pretty_log_
 
index_t level_
 
bool save_RVD_iter_
 
bool save_RVD_last_iter_
 
bool show_RVD_seed_
 
index_t current_iter_
 
bool newton_
 
bool verbose_
 
double epsilon_regularization_
 Add a regularization term to remove translational degree of freedom for the weights.
 
index_t nbZ_
 Number of empty cells in last iteration.
 
double g_norm_
 Norm of the gradient in last iteration.
 
double measure_of_smallest_cell_
 Measure of the smallest Laguerre cell.
 
bool w_did_not_change_
 True if w did not change, thus there is no need to recompute the power diagram.
 
double * Laguerre_centroids_
 If user-specified, then Laguerre centroids are output here.
 
double linsolve_epsilon_
 maximum value of \( \| Ax - b \| / \| b \| \)
 
index_t linsolve_maxiter_
 maximum number of iterations for linear solve
 
index_t linesearch_maxiter_
 maximum number of steplength divisions
 
index_t linesearch_init_iter_
 starting number of steplength divisions
 
OTLinearSolver linear_solver_
 one of OT_PRECG, OT_SUPERLU, OT_CHOLMOD.
 
const double * air_particles_
 if set, pointer to the air particles.
 
index_t nb_air_particles_
 if non-zero, number of air particles.
 
index_t air_particles_stride_
 Number of doubles between two consecutive air particles in air_particles_.
 
double air_fraction_
 The fraction of the total mass occupied by air.
 
bool clip_by_balls_
 Enabled if air fraction is specified without any air particles.
 
bool user_H_g_
 True if class is just used by user to compute Hessian and gradient instead of doing full computation.
 
NLMatrix user_H_
 User-defined Hessian matrix.
 
- Static Protected Attributes inherited from GEO::OptimalTransportMap
static OptimalTransportMapinstance_
 

Detailed Description

Computes semi-discrete optimal transport maps.

Computes an optimal transport map between two distributions in 2D. The first distribution is represented by a 2D triangulated mesh. The second distribution is a sum of Diracs with 2D coordinates. The algorithm is described in the following references:

  • 3D algorithm: http://arxiv.org/abs/1409.1279
  • Earlier 2D version by Quentin M\'erigot: Q. Merigot. A multiscale approach to optimal transport. Computer Graphics Forum 30 (5) 1583–1592, 2011 (Proc SGP 2011).
  • Earlier article on OT and power diagrams: F. Aurenhammer, F. Hoffmann, and B. Aronov. Minkowski-type theorems and least-squares clustering. Algorithmica, 20:61-76, 1998.

Definition at line 110 of file optimal_transport_2d.h.

Constructor & Destructor Documentation

◆ OptimalTransportMap2d()

GEO::OptimalTransportMap2d::OptimalTransportMap2d ( Mesh mesh,
const std::string &  delaunay = "BPOW2d",
bool  BRIO = false 
)

OptimalTransportMap2d constructor.

Parameters
[in]meshthe source distribution, represented as a 2d mesh. It can be also a 3D mesh with the Z coordinate set to 0.
[in]delaunayfactory name of the Delaunay triangulation.
[in]BRIOtrue if vertices are already ordered using BRIO

Member Function Documentation

◆ call_callback_on_RVD()

void GEO::OptimalTransportMap2d::call_callback_on_RVD ( )
overrideprotectedvirtual

Calls the callback for each intersection between a Laguerre cell and a simplex of the background mesh.

Implements GEO::OptimalTransportMap.

◆ compute_Laguerre_centroids()

void GEO::OptimalTransportMap2d::compute_Laguerre_centroids ( double *  centroids)
overridevirtual

Computes the centroids of the Laguerre cells.

Parameters
[out]centroidsa pointer to the dimension()*nb_points coordinates of the centroids.

Implements GEO::OptimalTransportMap.

◆ get_RVD()

void GEO::OptimalTransportMap2d::get_RVD ( Mesh M)
overridevirtual

Computes a mesh with the restricted Voronoi diagram.

Parameters
[out]Ma reference to the computed restricted Voronoi diagram.

Implements GEO::OptimalTransportMap.

◆ total_mesh_mass()

double GEO::OptimalTransportMap2d::total_mesh_mass ( ) const

Gets the total mass of the mesh.

Take the weights into account if they are present.

Returns
the total mass of the mesh.

Member Data Documentation

◆ clipped_

GEOGen::Polygon GEO::OptimalTransportMap2d::clipped_

Used by clipping operations.

Definition at line 161 of file optimal_transport_2d.h.

◆ work_

GEOGen::Polygon GEO::OptimalTransportMap2d::work_

Used by clipping operations.

Definition at line 157 of file optimal_transport_2d.h.


The documentation for this class was generated from the following file: