![]() |
SemiDiscreteOT 1.0
Semi-Discrete Optimal Transport Library
|
A class for refining the optimal transport potential using a softmax operation. More...
#include <SoftmaxRefinement.h>
Classes | |
struct | CopyData |
A struct to hold copy data for parallel assembly. More... | |
struct | ScratchData |
A struct to hold scratch data for parallel assembly. More... | |
Public Member Functions | |
void | set_distance_function (const std::function< double(const Point< spacedim > &, const Point< spacedim > &)> &dist) |
Sets the distance function to be used. | |
SoftmaxRefinement (MPI_Comm mpi_comm, const DoFHandler< dim, spacedim > &dof_handler, const Mapping< dim, spacedim > &mapping, const FiniteElement< dim, spacedim > &fe, const LinearAlgebra::distributed::Vector< double > &source_density, unsigned int quadrature_order, double distance_threshold, bool use_log_sum_exp_trick=true) | |
Constructor for the SoftmaxRefinement class. | |
Vector< double > | compute_refinement (const std::vector< Point< spacedim > > &target_points_fine, const Vector< double > &target_density_fine, const std::vector< Point< spacedim > > &target_points_coarse, const Vector< double > &target_density_coarse, const Vector< double > &potential_coarse, double regularization_param, int current_level, const std::vector< std::vector< std::vector< size_t > > > &child_indices) |
Computes the refined potential. | |
Private Types | |
using | IndexedPoint = std::pair< Point< spacedim >, std::size_t > |
using | RTreeParams = boost::geometry::index::rstar< 8 > |
using | RTree = boost::geometry::index::rtree< IndexedPoint, RTreeParams > |
Private Member Functions | |
void | setup_rtree () |
Sets up the R-tree. | |
std::vector< std::size_t > | find_nearest_target_points (const Point< spacedim > &query_point) const |
Finds the nearest target points to a query point. | |
void | local_assemble (const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, ScratchData &scratch_data, CopyData ©_data) |
Assembles the local contributions to the refined potential. | |
bool | is_simplex_element () const |
Checks if the finite element is a simplex element. | |
Private Attributes | |
MPI_Comm | mpi_communicator |
The MPI communicator. | |
const unsigned int | n_mpi_processes |
The number of MPI processes. | |
const unsigned int | this_mpi_process |
The rank of the current MPI process. | |
ConditionalOStream | pcout |
A conditional output stream for parallel printing. | |
std::function< double(const Point< spacedim > &, const Point< spacedim > &)> | distance_function |
The distance function. | |
const DoFHandler< dim, spacedim > & | dof_handler |
The DoF handler for the source mesh. | |
const Mapping< dim, spacedim > & | mapping |
The mapping for the source mesh. | |
const FiniteElement< dim, spacedim > & | fe |
The finite element for the source mesh. | |
const LinearAlgebra::distributed::Vector< double > & | source_density |
The density of the source measure. | |
const unsigned int | quadrature_order |
The order of the quadrature rule to use for integration. | |
RTree | target_points_rtree |
An R-tree for fast spatial queries on the target points. | |
double | current_lambda {0.0} |
The current regularization parameter. | |
const double | current_distance_threshold |
The current distance threshold. | |
const bool | use_log_sum_exp_trick |
Whether to use the log-sum-exp trick. | |
const std::vector< Point< spacedim > > * | current_target_points_fine {nullptr} |
The current fine target points. | |
const Vector< double > * | current_target_density_fine {nullptr} |
The current fine target density. | |
const std::vector< Point< spacedim > > * | current_target_points_coarse {nullptr} |
The current coarse target points. | |
const Vector< double > * | current_target_density_coarse {nullptr} |
The current coarse target density. | |
const Vector< double > * | current_potential_coarse {nullptr} |
The current coarse potential. | |
const std::vector< std::vector< std::vector< size_t > > > * | current_child_indices {nullptr} |
The current child indices. | |
int | current_level {0} |
The current level. | |
A class for refining the optimal transport potential using a softmax operation.
This class implements a method for refining the optimal transport potential from a coarser level to a finer level in a multilevel scheme. It uses a softmax-like operation to distribute the potential from the coarse points to their children in the finer level.
dim | The dimension of the source mesh. |
spacedim | The dimension of the space the mesh is embedded in. |
Definition at line 43 of file SoftmaxRefinement.h.
|
private |
Definition at line 147 of file SoftmaxRefinement.h.
|
private |
Definition at line 148 of file SoftmaxRefinement.h.
|
private |
Definition at line 149 of file SoftmaxRefinement.h.
SoftmaxRefinement< dim, spacedim >::SoftmaxRefinement | ( | MPI_Comm | mpi_comm, |
const DoFHandler< dim, spacedim > & | dof_handler, | ||
const Mapping< dim, spacedim > & | mapping, | ||
const FiniteElement< dim, spacedim > & | fe, | ||
const LinearAlgebra::distributed::Vector< double > & | source_density, | ||
unsigned int | quadrature_order, | ||
double | distance_threshold, | ||
bool | use_log_sum_exp_trick = true |
||
) |
Constructor for the SoftmaxRefinement class.
mpi_comm | The MPI communicator. |
dof_handler | The DoF handler for the source mesh. |
mapping | The mapping for the source mesh. |
fe | The finite element for the source mesh. |
source_density | The density of the source measure. |
quadrature_order | The order of the quadrature rule to use for integration. |
distance_threshold | The distance threshold for the R-tree search. |
use_log_sum_exp_trick | Whether to use the log-sum-exp trick for numerical stability. |
Definition at line 5 of file SoftmaxRefinement.cc.
|
inline |
Sets the distance function to be used.
dist | The distance function. |
Definition at line 49 of file SoftmaxRefinement.h.
Vector< double > SoftmaxRefinement< dim, spacedim >::compute_refinement | ( | const std::vector< Point< spacedim > > & | target_points_fine, |
const Vector< double > & | target_density_fine, | ||
const std::vector< Point< spacedim > > & | target_points_coarse, | ||
const Vector< double > & | target_density_coarse, | ||
const Vector< double > & | potential_coarse, | ||
double | regularization_param, | ||
int | current_level, | ||
const std::vector< std::vector< std::vector< size_t > > > & | child_indices | ||
) |
Computes the refined potential.
target_points_fine | The points of the fine target measure. |
target_density_fine | The weights of the fine target measure. |
target_points_coarse | The points of the coarse target measure. |
target_density_coarse | The weights of the coarse target measure. |
potential_coarse | The potential at the coarse level. |
regularization_param | The regularization parameter. |
current_level | The current level in the multilevel hierarchy. |
child_indices | The indices of the children of each coarse point. |
Definition at line 201 of file SoftmaxRefinement.cc.
|
private |
Sets up the R-tree.
Definition at line 30 of file SoftmaxRefinement.cc.
|
private |
Finds the nearest target points to a query point.
query_point | The query point. |
Definition at line 44 of file SoftmaxRefinement.cc.
|
private |
Assembles the local contributions to the refined potential.
cell | The current cell. |
scratch_data | The scratch data. |
copy_data | The copy data. |
Definition at line 62 of file SoftmaxRefinement.cc.
|
inlineprivate |
Checks if the finite element is a simplex element.
Definition at line 189 of file SoftmaxRefinement.h.
|
private |
The MPI communicator.
Definition at line 131 of file SoftmaxRefinement.h.
|
private |
The number of MPI processes.
Definition at line 132 of file SoftmaxRefinement.h.
|
private |
The rank of the current MPI process.
Definition at line 133 of file SoftmaxRefinement.h.
|
private |
A conditional output stream for parallel printing.
Definition at line 134 of file SoftmaxRefinement.h.
|
private |
The distance function.
Definition at line 137 of file SoftmaxRefinement.h.
|
private |
The DoF handler for the source mesh.
Definition at line 140 of file SoftmaxRefinement.h.
|
private |
The mapping for the source mesh.
Definition at line 141 of file SoftmaxRefinement.h.
|
private |
The finite element for the source mesh.
Definition at line 142 of file SoftmaxRefinement.h.
|
private |
The density of the source measure.
Definition at line 143 of file SoftmaxRefinement.h.
|
private |
The order of the quadrature rule to use for integration.
Definition at line 144 of file SoftmaxRefinement.h.
|
private |
An R-tree for fast spatial queries on the target points.
Definition at line 150 of file SoftmaxRefinement.h.
|
private |
The current regularization parameter.
Definition at line 153 of file SoftmaxRefinement.h.
|
private |
The current distance threshold.
Definition at line 154 of file SoftmaxRefinement.h.
|
private |
Whether to use the log-sum-exp trick.
Definition at line 155 of file SoftmaxRefinement.h.
|
private |
The current fine target points.
Definition at line 156 of file SoftmaxRefinement.h.
|
private |
The current fine target density.
Definition at line 157 of file SoftmaxRefinement.h.
|
private |
The current coarse target points.
Definition at line 158 of file SoftmaxRefinement.h.
|
private |
The current coarse target density.
Definition at line 159 of file SoftmaxRefinement.h.
|
private |
The current coarse potential.
Definition at line 160 of file SoftmaxRefinement.h.
|
private |
The current child indices.
Definition at line 161 of file SoftmaxRefinement.h.
|
private |
The current level.
Definition at line 162 of file SoftmaxRefinement.h.