![]() |
SemiDiscreteOT 1.0
Semi-Discrete Optimal Transport Library
|
Collection of utility functions for file I/O, mesh handling and data management. More...
Functions | |
std::string | to_scientific_string (double value, int precision=6) |
Convert a double to a string in scientific notation. | |
template<typename VectorContainer > | |
void | write_vector (const VectorContainer &points, const std::string &filepath, const std::string &fileMode="txt") |
Write a vector container to a file in binary or text format. | |
template<typename VectorContainer > | |
bool | read_vector (VectorContainer &points, const std::string &filepath, const std::string &fileMode="") |
Read a vector container from a file in binary or text format. | |
std::vector< std::string > | select_folder (const std::string &base_dir="output", bool allow_all=true) |
List available epsilon folders and allow user selection. | |
std::vector< std::pair< std::string, std::string > > | get_target_hierarchy_files (const std::string &dir) |
Get target point cloud hierarchy files. | |
template<int dim, int spacedim> | |
bool | load_hierarchy_data (const std::string &hierarchy_dir, std::vector< std::vector< std::vector< size_t > > > &child_indices, int specific_level, const MPI_Comm &mpi_communicator, dealii::ConditionalOStream &pcout) |
Load hierarchy data from files. | |
template<int dim, int spacedim = dim> | |
bool | write_mesh (const dealii::Triangulation< dim, spacedim > &mesh, const std::string &filepath, const std::vector< std::string > &formats, const std::vector< double > *cell_data=nullptr, const std::string &data_name="cell_data") |
Write mesh to file in specified formats with optional cell data. | |
template<int dim, int spacedim> | |
void | interpolate_non_conforming_nearest (const dealii::DoFHandler< dim, spacedim > &source_dh, const dealii::Vector< double > &source_field, const dealii::DoFHandler< dim, spacedim > &target_dh, dealii::LinearAlgebra::distributed::Vector< double > &target_field) |
Interpolate a field from a source mesh to a target mesh using nearest neighbor approach. | |
template<int dim, int spacedim> | |
void | interpolate_non_conforming (const dealii::DoFHandler< dim, spacedim > &source_dh, const dealii::Vector< double > &source_field, const dealii::DoFHandler< dim, spacedim > &target_dh, dealii::LinearAlgebra::distributed::Vector< double > &target_field) |
Interpolate a field from a source mesh to a target mesh. | |
template<int dim, int spacedim = dim> | |
std::unique_ptr< dealii::FiniteElement< dim, spacedim > > | create_fe_for_mesh (const dealii::Triangulation< dim, spacedim > &triangulation, const unsigned int degree=1) |
Detect cell types in a triangulation and create appropriate finite element. | |
template<int dim, int spacedim> | |
std::pair< std::unique_ptr< dealii::FiniteElement< dim, spacedim > >, std::unique_ptr< dealii::Mapping< dim, spacedim > > > | create_fe_and_mapping_for_mesh (const dealii::Triangulation< dim, spacedim > &triangulation, const unsigned int fe_degree=1, const unsigned int mapping_degree=1) |
Create appropriate finite element and mapping for a triangulation. | |
template<int dim, int spacedim = dim> | |
std::unique_ptr< dealii::Quadrature< dim > > | create_quadrature_for_mesh (const dealii::Triangulation< dim, spacedim > &triangulation, const unsigned int order=2) |
Create appropriate quadrature for a triangulation based on cell types. | |
template<int spacedim> | |
bool | write_points_with_density_vtk (const std::vector< dealii::Point< spacedim > > &points, const std::vector< double > &density, const std::string &filename, const std::string &description="Points with density values", const std::string &density_name="density") |
Write points with associated density to VTK file. | |
template<int spacedim> | |
bool | write_points_with_displacement_vtk (const std::vector< dealii::Point< spacedim > > &source_points, const std::vector< dealii::Point< spacedim > > &mapped_points, const std::vector< double > &density, const std::string &filename, const std::string &description="Points with displacement vectors", const std::string &density_name="density") |
Write points with displacement vectors and density to VTK file. | |
template<int dim, int spacedim = dim> | |
bool | read_vtk_field (const std::string &filename, dealii::DoFHandler< dim, spacedim > &vtk_dof_handler, dealii::Vector< double > &vtk_field, dealii::Triangulation< dim, spacedim > &vtk_tria, const MPI_Comm &mpi_communicator, dealii::ConditionalOStream &pcout, bool broadcast_field=false, const std::string &field_name="normalized_density") |
Read a scalar field from a VTK file. | |
Collection of utility functions for file I/O, mesh handling and data management.
|
inline |
void Utils::write_vector | ( | const VectorContainer & | points, |
const std::string & | filepath, | ||
const std::string & | fileMode = "txt" |
||
) |
Write a vector container to a file in binary or text format.
VectorContainer | Type of vector container |
points | Vector container to write |
filepath | Path to output file (without extension) |
fileMode | Output format ("txt" or "bin") |
Definition at line 56 of file utils.h.
bool Utils::read_vector | ( | VectorContainer & | points, |
const std::string & | filepath, | ||
const std::string & | fileMode = "" |
||
) |
Read a vector container from a file in binary or text format.
VectorContainer | Type of vector container |
points | Vector container to store read data |
filepath | Path to input file (with or without extension) |
fileMode | Input format ("txt" or "bin"), if empty will be inferred from filepath |
Definition at line 108 of file utils.h.
|
inline |
|
inline |
bool Utils::load_hierarchy_data | ( | const std::string & | hierarchy_dir, |
std::vector< std::vector< std::vector< size_t > > > & | child_indices, | ||
int | specific_level, | ||
const MPI_Comm & | mpi_communicator, | ||
dealii::ConditionalOStream & | pcout | ||
) |
Load hierarchy data from files.
hierarchy_dir | Directory containing hierarchy files |
child_indices | Vector to store parent-child relationships |
specific_level | Level to load (-1 for all levels) |
mpi_communicator | MPI communicator |
pcout | Parallel console output |
bool Utils::write_mesh | ( | const dealii::Triangulation< dim, spacedim > & | mesh, |
const std::string & | filepath, | ||
const std::vector< std::string > & | formats, | ||
const std::vector< double > * | cell_data = nullptr , |
||
const std::string & | data_name = "cell_data" |
||
) |
Write mesh to file in specified formats with optional cell data.
dim | Dimension of the mesh |
spacedim | Ambient space dimension |
mesh | Triangulation to write |
filepath | Base path for output files (without extension) |
formats | Vector of output formats ("vtk", "msh", "vtu") |
cell_data | Optional vector of cell data to include |
data_name | Name for the cell data field |
Definition at line 410 of file utils.h.
void Utils::interpolate_non_conforming_nearest | ( | const dealii::DoFHandler< dim, spacedim > & | source_dh, |
const dealii::Vector< double > & | source_field, | ||
const dealii::DoFHandler< dim, spacedim > & | target_dh, | ||
dealii::LinearAlgebra::distributed::Vector< double > & | target_field | ||
) |
Interpolate a field from a source mesh to a target mesh using nearest neighbor approach.
dim | Dimension of the mesh |
spacedim | Spatial dimension |
source_dh | Source DoFHandler |
source_field | Source field values |
target_dh | Target DoFHandler |
target_field | Target field values (output) |
Definition at line 504 of file utils.h.
void Utils::interpolate_non_conforming | ( | const dealii::DoFHandler< dim, spacedim > & | source_dh, |
const dealii::Vector< double > & | source_field, | ||
const dealii::DoFHandler< dim, spacedim > & | target_dh, | ||
dealii::LinearAlgebra::distributed::Vector< double > & | target_field | ||
) |
Interpolate a field from a source mesh to a target mesh.
dim | Dimension of the mesh |
spacedim | Spatial dimension |
source_dh | Source DoFHandler |
source_field | Source field values |
target_dh | Target DoFHandler |
target_field | Target field values (output) |
std::unique_ptr< dealii::FiniteElement< dim, spacedim > > Utils::create_fe_for_mesh | ( | const dealii::Triangulation< dim, spacedim > & | triangulation, |
const unsigned int | degree = 1 |
||
) |
Detect cell types in a triangulation and create appropriate finite element.
dim | Dimension of the mesh |
triangulation | Input triangulation to analyze |
degree | Polynomial degree for the finite element (defaults to 1) |
std::runtime_error | if cell type cannot be determined |
Definition at line 743 of file utils.h.
std::pair< std::unique_ptr< dealii::FiniteElement< dim, spacedim > >, std::unique_ptr< dealii::Mapping< dim, spacedim > > > Utils::create_fe_and_mapping_for_mesh | ( | const dealii::Triangulation< dim, spacedim > & | triangulation, |
const unsigned int | fe_degree = 1 , |
||
const unsigned int | mapping_degree = 1 |
||
) |
Create appropriate finite element and mapping for a triangulation.
dim | Dimension of the mesh |
triangulation | Input triangulation to analyze |
fe_degree | Polynomial degree for the finite element (defaults to 1) |
mapping_degree | Polynomial degree for the mapping (defaults to 1) |
std::runtime_error | if cell type cannot be determined |
std::unique_ptr< dealii::Quadrature< dim > > Utils::create_quadrature_for_mesh | ( | const dealii::Triangulation< dim, spacedim > & | triangulation, |
const unsigned int | order = 2 |
||
) |
Create appropriate quadrature for a triangulation based on cell types.
dim | Dimension of the mesh |
triangulation | Input triangulation to analyze |
order | Quadrature order (defaults to 2) |
std::runtime_error | if cell type cannot be determined |
bool Utils::write_points_with_density_vtk | ( | const std::vector< dealii::Point< spacedim > > & | points, |
const std::vector< double > & | density, | ||
const std::string & | filename, | ||
const std::string & | description = "Points with density values" , |
||
const std::string & | density_name = "density" |
||
) |
Write points with associated density to VTK file.
spacedim | Spatial dimension |
points | Vector of points to write |
density | Vector of density values associated with points |
filename | Output VTK filename |
description | Description for the VTK file header |
density_name | Name for the density scalar field |
bool Utils::write_points_with_displacement_vtk | ( | const std::vector< dealii::Point< spacedim > > & | source_points, |
const std::vector< dealii::Point< spacedim > > & | mapped_points, | ||
const std::vector< double > & | density, | ||
const std::string & | filename, | ||
const std::string & | description = "Points with displacement vectors" , |
||
const std::string & | density_name = "density" |
||
) |
Write points with displacement vectors and density to VTK file.
spacedim | Spatial dimension |
source_points | Vector of source points |
mapped_points | Vector of mapped points (targets) |
density | Vector of density values associated with points |
filename | Output VTK filename |
description | Description for the VTK file header |
density_name | Name for the density scalar field |
bool Utils::read_vtk_field | ( | const std::string & | filename, |
dealii::DoFHandler< dim, spacedim > & | vtk_dof_handler, | ||
dealii::Vector< double > & | vtk_field, | ||
dealii::Triangulation< dim, spacedim > & | vtk_tria, | ||
const MPI_Comm & | mpi_communicator, | ||
dealii::ConditionalOStream & | pcout, | ||
bool | broadcast_field = false , |
||
const std::string & | field_name = "normalized_density" |
||
) |
Read a scalar field from a VTK file.
filename | Path to the VTK file |
vtk_dof_handler | DoF handler for the VTK mesh |
vtk_field | Vector for the VTK field data |
vtk_tria | Triangulation for the VTK mesh |
mpi_communicator | MPI communicator |
pcout | Parallel output stream |
broadcast_field | Whether to broadcast the field to all processes (true for source, false for target) |
field_name | Name of the field to read from the VTK file |
Definition at line 1039 of file utils.h.