FEAT 3
Finite Element Analysis Toolbox
Loading...
Searching...
No Matches
FEAT::Geometry::VoxelMap Class Reference

VoxelMap class. More...

#include <voxel_map.hpp>

Classes

struct  FileHeader
 header structure of voxel map file header More...
 
class  ReadResult
 helper class for read function result More...
 
class  WriteResult
 helper class for write function result More...
 

Public Types

typedef std::int64_t i64
 signed 64-bit integer type More...
 
typedef std::uint64_t u64
 unsigned 64-bit integer type More...
 

Public Member Functions

 VoxelMap ()
 standard constructor More...
 
 VoxelMap (const VoxelMap &)=delete
 no copies, no problems
 
 VoxelMap (VoxelMap &&other)
 move constructor More...
 
virtual ~VoxelMap ()
 virtual destructor More...
 
template<typename Coord_ , int dim_>
bool check_box (const Tiny::Vector< Coord_, dim_ > &box_min, const Tiny::Vector< Coord_, dim_ > &box_max) const
 Checks the voxel map entries for a given box of voxels. More...
 
bool check_point (u64 ix, u64 iy=0u, u64 iz=0u) const
 Returns the value of a specific voxel point in the voxel map. More...
 
template<typename Coord_ , int dim_>
bool check_point_nearest (const Tiny::Vector< Coord_, dim_ > &point) const
 Checks the voxel map entry for a given point. More...
 
template<typename Coord_ , int dim_>
void compute_map (const Dist::Comm &comm, VoxelMasker< Coord_, dim_ > &masker, bool gather_to_all=true)
 Creates the voxel map based on a given masker object by utilizing MPI parallelization. More...
 
template<typename MeshType_ >
void compute_map_from_chart (const Dist::Comm &comm, const Geometry::Atlas::ChartBase< MeshType_ > &chart, bool invert, bool gather_to_all=true)
 Creates the voxel map based on a chart by utilizing MPI parallelization. More...
 
void compute_map_from_formula_2d (const Dist::Comm &comm, const String &formula, bool gather_to_all=true)
 Creates the voxel map based on a 2D formula by utilizing MPI parallelization. More...
 
void compute_map_from_formula_3d (const Dist::Comm &comm, const String &formula, bool gather_to_all=true)
 Creates the voxel map based on a 3D formula by utilizing MPI parallelization. More...
 
template<typename Lambda_ >
void compute_map_from_lambda_2d (const Dist::Comm &comm, Lambda_ &&lambda, bool gather_to_all=true)
 Creates the voxel map based on a 2D lambda expression by utilizing MPI parallelization. More...
 
template<typename Lambda_ >
void compute_map_from_lambda_3d (const Dist::Comm &comm, Lambda_ &&lambda, bool gather_to_all=true)
 Creates the voxel map based on a 3D lambda expression by utilizing MPI parallelization. More...
 
void compute_map_from_off_3d (const Dist::Comm &comm, const String &filename, bool invert, bool gather_to_all=true)
 Creates the voxel map based on a 3D OFF model handled by CGAL by utilizing MPI parallelization. More...
 
void compute_map_from_off_3d (const Dist::Comm &comm, std::istream &is, bool invert, bool gather_to_all=true)
 Creates the voxel map based on a 3D OFF model handled by CGAL by utilizing MPI parallelization. More...
 
void export_plane_to_bmp (const String &filename, Real z_coord) const
 Exports a single plane of the voxel map to a monochrome BMP image file. More...
 
void export_plane_to_bmp (std::ostream &os, Real z_coord) const
 Exports a single plane of the voxel map to a monochrome BMP image file. More...
 
void export_to_bmp (const String &filename_prefix) const
 Exports the voxel map to a sequence of monochrome BMP image files. More...
 
Real get_bounding_box_max (int dim) const
 Returns the maximum bounding box coordinate for a given dimensions. More...
 
Real get_bounding_box_min (int dim) const
 Returns the minimum bounding box coordinate for a given dimensions. More...
 
Real get_bounding_box_volume () const
 Returns the volume of the 1D/2D/3D voxel map bounding box. More...
 
Real get_domain_coverage () const
 Returns the domain coverage of the voxel map. More...
 
Real get_domain_volume () const
 Returns the covered domain volume. More...
 
const std::vector< char > & get_map () const
 Returns a reference to the voxel map. More...
 
Index get_num_points (int dim) const
 Returns the number of points of the voxel map in each dimension. More...
 
Index get_num_voxels () const
 Returns the total number of voxels in the voxel map. More...
 
bool get_out_of_bounds_value () const
 Returns the out-of-bounds value for the voxel map. More...
 
u64 get_stride_line () const
 Returns the line stride, i.e. the size of a single X-line in bytes. More...
 
u64 get_stride_plane () const
 Returns the plane stride, i.e. the size of a single XY-plane in bytes. More...
 
u64 get_stride_volume () const
 Returns the volume stride, i.e. the size of a single XYZ-volume in bytes. More...
 
VoxelMapoperator= (const VoxelMap &)=delete
 
VoxelMapoperator= (VoxelMap &&other)
 move assign operator More...
 
ReadResult read (const Dist::Comm &comm, const String &filename)
 Reads a voxel map from a file. More...
 
ReadResult read (std::istream &is, const String &filename="")
 Reads a voxel map from a binary input stream. More...
 
void render_plane_to_bmp (const String &filename, Index width, Index height, Real z_min, Real z_max) const
 Renders a plane range of voxel map to a single gray-scale BMP image file. More...
 
void render_plane_to_bmp (std::ostream &os, Index width, Index height, Real z_min, Real z_max) const
 Renders a plane range of voxel map to a single gray-scale BMP image file stream. More...
 
void render_to_bmp (const String &filename_prefix, Index width, Index height, Index depth) const
 Renders the voxel map to a sequence of gray-scale BMP image files. More...
 
template<typename Coord_ , int dim_>
Real sample_box (const Tiny::Vector< Coord_, dim_ > &box_min, const Tiny::Vector< Coord_, dim_ > &box_max) const
 Samples the voxel map entries for a given box of voxels. More...
 
template<typename Coord_ , int dim_>
Real sample_point (const Tiny::Vector< Coord_, dim_ > &point) const
 Samples the voxel map entry for a given point by multi-linear interpolation. More...
 
void set_bounding_box_2d (Real x_min, Real x_max, Real y_min, Real y_max)
 Sets the bounding box for a 2D voxel map. More...
 
void set_bounding_box_3d (Real x_min, Real x_max, Real y_min, Real y_max, Real z_min, Real z_max)
 Sets the bounding box for a 3D voxel map. More...
 
void set_num_points (Index num_x, Index num_y, Index num_z=0u)
 Sets the number of points of the voxel map in each dimension. More...
 
void set_out_of_bounds_value (bool value)
 Sets the out-of-bounds value for the voxel map. More...
 
void set_resolution (Real max_res)
 Sets the resolution for the voxel map, i.e. the maximum voxel size in each dimension. More...
 
WriteResult write (const String &filename, const u64 compress_block_size=128u, const int compress_level=9) const
 Writes the voxel map into a file that can later be read in. More...
 
WriteResult write (std::ostream &os, const u64 compress_block_size=128u, const int compress_level=9) const
 Writes the voxel map into a binary output stream. More...
 

Static Public Member Functions

static u64 calc_line_stride (u64 num_x)
 ensure that we're not compiling on some freak architecture... More...
 

Static Public Attributes

static constexpr u64 header_size = 136ull
 size of voxel map header in bytes More...
 
static constexpr u64 magic_number = 0x70614D6C65786F56ull
 voxel map magic number More...
 
static constexpr i64 unit_size = 1'000'000'000
 size of a real unit (meter) in our internal units (nanometers) More...
 

Protected Member Functions

Real _calc_sample_rate (i64 xyz, i64 idx, std::size_t sdim) const
 Helper function for _sample_voxel_map_Xd: compute the sample rate in a given dimension. More...
 
bool _check_box (const std::array< u64, 1 > &box_min, const std::array< u64, 1 > &box_max) const
 Checks whether a 1D box contains at least one active voxel. More...
 
bool _check_box (const std::array< u64, 2 > &box_min, const std::array< u64, 2 > &box_max) const
 Checks whether a 2D box contains at least one active voxel. More...
 
bool _check_box (const std::array< u64, 3 > &box_min, const std::array< u64, 3 > &box_max) const
 Checks whether a 3D box contains at least one active voxel. More...
 
bool _check_point (i64 xidx, i64 yidx, i64 zidx) const
 Helper function for _sample_point: checks a single voxel point value. More...
 
bool _check_point_nearest (const std::array< i64, 1 > &p) const
 Checks the nearest voxel map entry for a given 1D point. More...
 
bool _check_point_nearest (const std::array< i64, 2 > &p) const
 Checks the nearest voxel map entry for a given 2D point. More...
 
bool _check_point_nearest (const std::array< i64, 3 > &p) const
 Checks the nearest voxel map entry for a given 3D point. More...
 
void _compress_voxel_map_line (const std::vector< int > &mask, const u64 line)
 Compresses a voxel map line into the voxel map. More...
 
u64 _compute_domain_coverage () const
 Computes the domain coverage of the voxel map. More...
 
template<typename CoordType , int dim_>
void _compute_voxel_map (const Dist::Comm &comm, VoxelMasker< CoordType, dim_ > &masker, bool gather_to_all)
 Computes the voxel map for the domain. More...
 
template<typename CoordType_ >
void _compute_voxel_map_lines (VoxelMasker< CoordType_, 1 > &masker, u64 beg, u64 end, u64 offset)
 Computes a range of voxel map lines for a 1D domain. More...
 
template<typename CoordType_ >
void _compute_voxel_map_lines (VoxelMasker< CoordType_, 2 > &masker, u64 beg, u64 end, u64 offset)
 Computes a range of voxel map lines for a 2D domain. More...
 
template<typename CoordType_ >
void _compute_voxel_map_lines (VoxelMasker< CoordType_, 3 > &masker, u64 beg, u64 end, u64 offset)
 Computes a range of voxel map lines for a 3D domain. More...
 
void _export_plane_to_bmp (std::ostream &os, u64 plane) const
 Exports a single plane of the voxel map to a monochrome BMP image file. More...
 
i64 _map_coord_idx_lower (i64 xyz, std::size_t sdim) const
 Maps a unit coordinate to a X-/Y-/Z-index of the lower voxel in the voxel map. More...
 
u64 _map_coord_idx_nearest (i64 xyz, std::size_t sdim) const
 Maps a unit coordinate to a X-/Y-/Z-index of the closest voxel in the voxel map. More...
 
i64 _map_coord_idx_upper (i64 xyz, std::size_t sdim) const
 Maps a unit coordinate to a X-/Y-/Z-index of the upper voxel in the voxel map. More...
 
void _render_plane_to_bmp (std::ostream &os, Index width, Index height, i64 box_min_z, i64 box_max_z) const
 Renders a plane range of voxel map to a single gray-scale BMP image file. More...
 
Real _sample_box (const std::array< i64, 1 > &box_min, const std::array< i64, 1 > &box_max) const
 Samples a 1D box of the voxel map. More...
 
Real _sample_box (const std::array< i64, 2 > &box_min, const std::array< i64, 2 > &box_max) const
 Samples a 2D box of the voxel map. More...
 
Real _sample_box (const std::array< i64, 3 > &box_min, const std::array< i64, 3 > &box_max) const
 Samples a 3D box of the voxel map. More...
 
Real _sample_box_1d (i64 px0, i64 px1, i64 xidx0, i64 xidx1, i64 yidx=i64(0), i64 zidx=i64(0)) const
 Helper function for _sample_box: samples a 1D edge in X line. More...
 
Real _sample_box_2d (i64 px0, i64 px1, i64 py0, i64 py1, i64 xidx0, i64 xidx1, i64 yidx0, i64 yidx1, i64 zidx=i64(0)) const
 Helper function for _sample_box: samples a 2D rectangle in XY plane. More...
 
Real _sample_box_3d (i64 px0, i64 px1, i64 py0, i64 py1, i64 pz0, i64 pz1, i64 xidx0, i64 xidx1, i64 yidx0, i64 yidx1, i64 zidx0, i64 zidx1) const
 Helper function for _sample_box: samples a 3D cuboid in XYZ volume. More...
 
Real _sample_point (const std::array< i64, 1 > &p) const
 Samples the voxel map for a given 1D point by multi-linear interpolation. More...
 
Real _sample_point (const std::array< i64, 2 > &p) const
 Samples the voxel map for a given 2D point by multi-linear interpolation. More...
 
Real _sample_point (const std::array< i64, 3 > &p) const
 Samples the voxel map for a given 3D point by multi-linear interpolation. More...
 
Real _sample_point_1d_x (i64 px, i64 xidx, i64 yidx=i64(0), i64 zidx=i64(0)) const
 Helper function for _sample_point: samples a 1D edge in X-parallel line. More...
 
Real _sample_point_1d_y (i64 py, i64 xidx, i64 yidx, i64 zidx=i64(0)) const
 Helper function for _sample_point: samples a 1D edge in Y-parallel line. More...
 
Real _sample_point_1d_z (i64 pz, i64 xidx, i64 yidx, i64 zidx) const
 Helper function for _sample_point: samples a 1D edge in Z-parallel line. More...
 
Real _sample_point_2d (i64 px, i64 py, i64 xidx, i64 yidx, i64 zidx=i64(0)) const
 Helper function for _sample_point: samples a 2D square in XY-parallel plane. More...
 
Real _sample_point_3d (i64 px, i64 py, i64 pz, i64 xidx, i64 yidx, i64 zidx) const
 Helper function for _sample_point: samples a 3D cube in XYZ volume. More...
 

Protected Attributes

std::array< i64, 3u > _bbox_max
 
std::array< i64, 3u > _bbox_min
 bounding box dimensions of domain More...
 
u64 _coverage
 the coverage of the domain, relative to unit size More...
 
int _create_stage
 current creation stage More...
 
u64 _num_lines
 number of lines; 2D = num_points[2]; 3D: = num_points[2] * num_points[3] More...
 
u64 _num_planes
 number of planes; 2D = 1; 3D = num_points[3] More...
 
std::array< u64, 3u > _num_points
 number of points in each dimension More...
 
bool _out_of_bounds_value
 the out-of-bounds-value for the voxel map More...
 
u64 _stride_line
 stride of a single voxel line in X dimension More...
 
u64 _stride_plane
 stride of a single plane in XY dimensions More...
 
std::vector< char > _voxel_map
 the actual voxel map; its size may be larger than necessary due to padding More...
 

Detailed Description

VoxelMap class.

This class is responsible for creating an managing voxel maps, which are used by the Control::Domain::VoxelDomainControl class to determine which elements can be dropped during the mesh hierarchy creation process.

A voxel map is basically a 2D or 3D monochrome raster image representing a domain bounding box, where each voxel in the map is either 1 or 0, depending on whether the point at the voxel's position is considered to be inside or outside of the domain that is supposed to be discretized.

Voxel maps can be created from analytic expressions which are either given as compile-time lambda expressions or as runtime string expressions to be parsed or from a surface triangulation given as an OFF file, assuming that the required third-party libraries 'CGAL' and 'fparser' are present.

This class also provides the functionality to write a voxel map into a binary output file, optionally compressed by the ZLIB library, along with the functionality to read these binary files, of course. It is strongly recommended to write voxel map in compressed format, because the compression often shrinks the resulting file sizes by more than 95%.

To create a voxel map object, one has to perform the following steps (in this particular order):

  1. Specify the bounding box for the voxel map by calling one of the following functions:
  2. Specify the number of points for the voxel map in each dimension by calling one of the following functions:
    • set_num_points to explicitly specify the number of points for each dimension
    • set_resolution to compute the number of points for a maximum voxel resolution
  3. Compute the voxel map contents by calling one of the following functions:
  4. Optional: Specify the out-of-bounds value by calling set_out_of_bounds_value

Peter Zajac

Definition at line 404 of file voxel_map.hpp.

Member Typedef Documentation

◆ i64

typedef std::int64_t FEAT::Geometry::VoxelMap::i64

signed 64-bit integer type

Definition at line 408 of file voxel_map.hpp.

◆ u64

typedef std::uint64_t FEAT::Geometry::VoxelMap::u64

unsigned 64-bit integer type

Definition at line 411 of file voxel_map.hpp.

Constructor & Destructor Documentation

◆ VoxelMap() [1/2]

FEAT::Geometry::VoxelMap::VoxelMap ( )

standard constructor

Definition at line 23 of file voxel_map.cpp.

◆ VoxelMap() [2/2]

FEAT::Geometry::VoxelMap::VoxelMap ( VoxelMap &&  other)

move constructor

Definition at line 38 of file voxel_map.cpp.

◆ ~VoxelMap()

FEAT::Geometry::VoxelMap::~VoxelMap ( )
virtual

virtual destructor

Definition at line 75 of file voxel_map.cpp.

Member Function Documentation

◆ _calc_sample_rate()

Real FEAT::Geometry::VoxelMap::_calc_sample_rate ( i64  xyz,
i64  idx,
std::size_t  sdim 
) const
protected

Helper function for _sample_voxel_map_Xd: compute the sample rate in a given dimension.

The sample rate of a point xyz with respect to the voxel v(idx) is defined as being equal to 1 if xyz = V(idx) and to fall linearly to 0 for xyz = V(idx-1) and xyz = V(idx+1), which can be achieved by the formula

                   sr(xyz,idx) :=  1 - |xyz - V(idx)| / |V(idx+1) - V(idx)|

note that the denominator |V(idx+1) - V(idx)| := 1/h is identical for all idx

Parameters
[in]xyzThe X-/Y-/Z unit coordinate of the point for which the sample rate is to be computed
[in]idxThe voxel index for the coordinate xyz as returned by _map_coord_idx_lower
Returns
The sample rate of the point xyz w.r.t. the voxel V(idx)

Definition at line 601 of file voxel_map.cpp.

References _bbox_min, _num_points, and FEAT::Math::abs().

Referenced by _sample_box_1d(), _sample_box_2d(), _sample_box_3d(), _sample_point_1d_x(), _sample_point_1d_y(), _sample_point_1d_z(), _sample_point_2d(), and _sample_point_3d().

◆ _check_box() [1/3]

bool FEAT::Geometry::VoxelMap::_check_box ( const std::array< u64, 1 > &  box_min,
const std::array< u64, 1 > &  box_max 
) const
protected

Checks whether a 1D box contains at least one active voxel.

Parameters
[in]box_min,box_maxThe minimum and maximum X indices that make up the box to be checked
Returns
true, if at least one voxel in the given box is set to true, otherwise false.

Definition at line 659 of file voxel_map.cpp.

References _voxel_map.

Referenced by check_box().

◆ _check_box() [2/3]

bool FEAT::Geometry::VoxelMap::_check_box ( const std::array< u64, 2 > &  box_min,
const std::array< u64, 2 > &  box_max 
) const
protected

Checks whether a 2D box contains at least one active voxel.

Parameters
[in]box_min,box_maxThe minimum and maximum X/Y indices that make up the box to be checked
Returns
true, if at least one voxel in the given box is set to true, otherwise false.

Definition at line 669 of file voxel_map.cpp.

References _stride_line, and _voxel_map.

◆ _check_box() [3/3]

bool FEAT::Geometry::VoxelMap::_check_box ( const std::array< u64, 3 > &  box_min,
const std::array< u64, 3 > &  box_max 
) const
protected

Checks whether a 3D box contains at least one active voxel.

Parameters
[in]box_min,box_maxThe minimum and maximum X/Y/Z indices that make up the box to be checked
Returns
true, if at least one voxel in the given box is set to true, otherwise false.

Definition at line 682 of file voxel_map.cpp.

References _num_points, _stride_line, and _voxel_map.

◆ _check_point()

bool FEAT::Geometry::VoxelMap::_check_point ( i64  xidx,
i64  yidx,
i64  zidx 
) const
protected

Helper function for _sample_point: checks a single voxel point value.

Definition at line 620 of file voxel_map.cpp.

References _num_points, _out_of_bounds_value, _stride_line, _stride_plane, and _voxel_map.

Referenced by _sample_box_1d(), and check_point().

◆ _check_point_nearest() [1/3]

bool FEAT::Geometry::VoxelMap::_check_point_nearest ( const std::array< i64, 1 > &  p) const
protected

Checks the nearest voxel map entry for a given 1D point.

Parameters
[in]pThe unit coordinates of the point whose closest voxel in the voxel map is to be checked.
Returns
true, if the closest voxel to the point is inside the voxel map domain (i.e. the voxel map entry is 1), or false, if the closest voxel is outside of the voxel map domain.

Definition at line 631 of file voxel_map.cpp.

References _map_coord_idx_nearest(), _out_of_bounds_value, and _voxel_map.

Referenced by check_point_nearest().

◆ _check_point_nearest() [2/3]

bool FEAT::Geometry::VoxelMap::_check_point_nearest ( const std::array< i64, 2 > &  p) const
protected

Checks the nearest voxel map entry for a given 2D point.

Parameters
[in]pThe unit coordinates of the point whose closest voxel in the voxel map is to be checked.
Returns
true, if the closest voxel to the point is inside the voxel map domain (i.e. the voxel map entry is 1), or false, if the closest voxel is outside of the voxel map domain.

Definition at line 639 of file voxel_map.cpp.

References _map_coord_idx_nearest(), _out_of_bounds_value, _stride_line, and _voxel_map.

◆ _check_point_nearest() [3/3]

bool FEAT::Geometry::VoxelMap::_check_point_nearest ( const std::array< i64, 3 > &  p) const
protected

Checks the nearest voxel map entry for a given 3D point.

Parameters
[in]pThe unit coordinates of the point whose closest voxel in the voxel map is to be checked.
Returns
true, if the closest voxel to the point is inside the voxel map domain (i.e. the voxel map entry is 1), or false, if the closest voxel is outside of the voxel map domain.

Definition at line 648 of file voxel_map.cpp.

References _map_coord_idx_nearest(), _num_points, _out_of_bounds_value, _stride_line, and _voxel_map.

◆ _compress_voxel_map_line()

void FEAT::Geometry::VoxelMap::_compress_voxel_map_line ( const std::vector< int > &  mask,
const u64  line 
)
protected

Compresses a voxel map line into the voxel map.

Parameters
[in]maskThe voxel map line as computed by a VoxelMasker implementation
[in]lineThe index of the line that is to be compressed.

Definition at line 997 of file voxel_map.cpp.

References _num_lines, _stride_line, _voxel_map, and ASSERTM.

Referenced by _compute_voxel_map_lines().

◆ _compute_domain_coverage()

u64 FEAT::Geometry::VoxelMap::_compute_domain_coverage ( ) const
protected

Computes the domain coverage of the voxel map.

Definition at line 518 of file voxel_map.cpp.

References _voxel_map, get_num_voxels(), and unit_size.

Referenced by _compute_voxel_map().

◆ _compute_voxel_map()

template<typename CoordType , int dim_>
void FEAT::Geometry::VoxelMap::_compute_voxel_map ( const Dist::Comm comm,
VoxelMasker< CoordType, dim_ > &  masker,
bool  gather_to_all 
)
inlineprotected

Computes the voxel map for the domain.

In an MPI-parallel simulation, this function splits the computation of the mask across all processes in the communicator and gathers the mask on all processes afterwards. In the borderline case where there are more processes in the communicator than there are voxel map lines in the mask, this function creates a temporary sub-communicator with as many processes as there are voxel map lines to distribute the work over these processes and performs a broadcast on the original communicator after gathering the mask on rank 0.

Parameters
[in]commA transient reference to the communicator that is to be used to spread the workload onto.
[in]maskerA transient reference to the masker that is to be used to compute the mask. Must be the same on all MPI processes
[in]gather_to_allSpecifies whether the voxel map is to be gathered on all processes in the communicator. The only reason to set this to false is when the map is only required on rank 0 for file output.

Definition at line 1623 of file voxel_map.hpp.

References _compute_domain_coverage(), _compute_voxel_map_lines(), FEAT::Dist::Comm::allgather(), FEAT::Dist::Comm::bcast(), FEAT::Dist::Comm::comm_create_range_incl(), FEAT::Dist::Comm::gather(), FEAT::Dist::Comm::is_null(), FEAT::Math::min(), FEAT::Dist::Comm::rank(), and FEAT::Dist::Comm::size().

Referenced by compute_map(), compute_map_from_chart(), compute_map_from_formula_2d(), compute_map_from_formula_3d(), compute_map_from_lambda_2d(), compute_map_from_lambda_3d(), and compute_map_from_off_3d().

◆ _compute_voxel_map_lines() [1/3]

template<typename CoordType_ >
void FEAT::Geometry::VoxelMap::_compute_voxel_map_lines ( VoxelMasker< CoordType_, 1 > &  masker,
u64  beg,
u64  end,
u64  offset 
)
inlineprotected

Computes a range of voxel map lines for a 1D domain.

Parameters
[in]maskerA transient reference to the masker that is to be used to compute the mask
[in]begThe index of the first voxel map line that is to be computed. Must be equal to 0 in 1D.
[in]endThe index of the last voxel map line that is to be computed plus one. Must be equal to 1 in 1D.
[in]offsetThe offset of the first voxel map line that is to be computed; Must be equal to 0 in 1D.

Definition at line 1510 of file voxel_map.hpp.

References _bbox_min, _compress_voxel_map_line(), _num_points, FEAT::Geometry::VoxelMasker< CoordType_, dim_ >::mask_line(), unit_size, and XASSERTM.

Referenced by _compute_voxel_map().

◆ _compute_voxel_map_lines() [2/3]

template<typename CoordType_ >
void FEAT::Geometry::VoxelMap::_compute_voxel_map_lines ( VoxelMasker< CoordType_, 2 > &  masker,
u64  beg,
u64  end,
u64  offset 
)
inlineprotected

Computes a range of voxel map lines for a 2D domain.

Parameters
[in]maskerA transient reference to the masker that is to be used to compute the mask
[in]begThe index of the first voxel map line that is to be computed.
[in]endThe index of the last voxel map line that is to be computed plus one.
[in]offsetThe offset of the first voxel map line that is to be computed; this is typically equal to 0, unless we are only computing a sub-map that is going to be gathered on the root later on, in which case it may be set equal to beg.

Definition at line 1540 of file voxel_map.hpp.

References _bbox_min, _compress_voxel_map_line(), _num_points, FEAT::Geometry::VoxelMasker< CoordType_, dim_ >::mask_line(), unit_size, and XASSERTM.

◆ _compute_voxel_map_lines() [3/3]

template<typename CoordType_ >
void FEAT::Geometry::VoxelMap::_compute_voxel_map_lines ( VoxelMasker< CoordType_, 3 > &  masker,
u64  beg,
u64  end,
u64  offset 
)
inlineprotected

Computes a range of voxel map lines for a 3D domain.

Parameters
[in]maskerA transient reference to the masker that is to be used to compute the mask
[in]begThe index of the first voxel map line that is to be computed.
[in]endThe index of the last voxel map line that is to be computed plus one.
[in]offsetThe offset of the first voxel map line that is to be computed; this is typically equal to 0, unless we are only computing a sub-map that is going to be gathered on the root later on, in which case it may be set equal to beg.

Definition at line 1576 of file voxel_map.hpp.

References _bbox_min, _compress_voxel_map_line(), FEAT::Geometry::VoxelMasker< CoordType_, dim_ >::mask_line(), and unit_size.

◆ _export_plane_to_bmp()

void FEAT::Geometry::VoxelMap::_export_plane_to_bmp ( std::ostream &  os,
u64  plane 
) const
protected

Exports a single plane of the voxel map to a monochrome BMP image file.

This function will export an entire plane of the voxel map in full resolution as a monochrome bitmap files.

Parameters
[in]filenameThe filename for the BMP image file
[in]planeThe index of the plane that is to be exported.

Definition at line 1008 of file voxel_map.cpp.

References _num_points, _stride_line, _stride_plane, _voxel_map, and XASSERTM.

Referenced by export_plane_to_bmp(), and export_to_bmp().

◆ _map_coord_idx_lower()

i64 FEAT::Geometry::VoxelMap::_map_coord_idx_lower ( i64  xyz,
std::size_t  sdim 
) const
protected

Maps a unit coordinate to a X-/Y-/Z-index of the lower voxel in the voxel map.

Attention
This function intentionally does not check the point against the valid voxel map bounding box, i.e. the returned index may be outside of the interval {0, ..., num_points[sdim]-1}.
Parameters
[in]xyzThe X-/Y-/Z unit coordinate that is to be mapped
[in]sdimSpecifies which dimension is to be mapped: 0=X, 1=Y, 2=Z
Returns
The index of the lower voxel to the given unit coordinate, i.e. the voxel whose unit coordinate is less than or equal to the given point unit coordinate xyz.

Definition at line 565 of file voxel_map.cpp.

References _bbox_min, and _num_points.

Referenced by _sample_box(), _sample_point(), and check_box().

◆ _map_coord_idx_nearest()

u64 FEAT::Geometry::VoxelMap::_map_coord_idx_nearest ( i64  xyz,
std::size_t  sdim 
) const
protected

Maps a unit coordinate to a X-/Y-/Z-index of the closest voxel in the voxel map.

Parameters
[in]xyzThe X-/Y-/Z unit coordinate that is to be mapped
[in]sdimSpecifies which dimension is to be mapped: 0=X, 1=Y, 2=Z
Returns
The index of the closest voxel to the given unit coordinate or ~u64(0), if the unit coordinate was outside of the voxel map bounding box.

Definition at line 535 of file voxel_map.cpp.

References _bbox_min, and _num_points.

Referenced by _check_point_nearest(), and export_plane_to_bmp().

◆ _map_coord_idx_upper()

i64 FEAT::Geometry::VoxelMap::_map_coord_idx_upper ( i64  xyz,
std::size_t  sdim 
) const
protected

Maps a unit coordinate to a X-/Y-/Z-index of the upper voxel in the voxel map.

Attention
This function intentionally does not check the point against the valid voxel map bounding box, i.e. the returned index may be outside of the interval {0, ..., num_points[sdim]-1}.
Parameters
[in]xyzThe X-/Y-/Z unit coordinate that is to be mapped
[in]sdimSpecifies which dimension is to be mapped: 0=X, 1=Y, 2=Z
Returns
The index of the upper voxel to the given unit coordinate, i.e. the voxel whose unit coordinate is greater than or equal to the given point unit coordinate xyz.

Definition at line 580 of file voxel_map.cpp.

References _bbox_min, and _num_points.

Referenced by _sample_box(), and check_box().

◆ _render_plane_to_bmp()

void FEAT::Geometry::VoxelMap::_render_plane_to_bmp ( std::ostream &  os,
Index  width,
Index  height,
i64  box_min_z,
i64  box_max_z 
) const
protected

Renders a plane range of voxel map to a single gray-scale BMP image file.

This function will render a range of plane of the voxel map in reduced resolution as a single gray-scale bitmap file, where each pixel contains the average value of all voxels that have been reduced to it.

Parameters
[in]osThe binary output stream for the BMP image
[in]widthThe desired bitmap width; must be less than or equal to the number of points in X-dimension.
[in]heightThe desired bitmap height; must be less than or equal to the number of points in Y-dimension.
[in]box_min_z,box_max_zThe indices of the plane Z range that has to be sampled into a single output plane.

Definition at line 1087 of file voxel_map.cpp.

References _bbox_min, _sample_box(), and XASSERTM.

Referenced by render_plane_to_bmp(), and render_to_bmp().

◆ _sample_box() [1/3]

Real FEAT::Geometry::VoxelMap::_sample_box ( const std::array< i64, 1 > &  box_min,
const std::array< i64, 1 > &  box_max 
) const
protected

Samples a 1D box of the voxel map.

Parameters
[in]box_min,box_maxThe minimum and maximum X unit coordinates that make up the box to be sampled
Returns
The average voxel map value of the sample box

Definition at line 957 of file voxel_map.cpp.

References _map_coord_idx_lower(), _map_coord_idx_upper(), _num_points, and _sample_box_1d().

Referenced by _render_plane_to_bmp(), and sample_box().

◆ _sample_box() [2/3]

Real FEAT::Geometry::VoxelMap::_sample_box ( const std::array< i64, 2 > &  box_min,
const std::array< i64, 2 > &  box_max 
) const
protected

Samples a 2D box of the voxel map.

Parameters
[in]box_min,box_maxThe minimum and maximum X/Y unit coordinates that make up the box to be sampled
Returns
The average voxel map value of the sample box

Definition at line 966 of file voxel_map.cpp.

References _map_coord_idx_lower(), _map_coord_idx_upper(), _num_points, _sample_box_1d(), and _sample_box_2d().

◆ _sample_box() [3/3]

Real FEAT::Geometry::VoxelMap::_sample_box ( const std::array< i64, 3 > &  box_min,
const std::array< i64, 3 > &  box_max 
) const
protected

Samples a 3D box of the voxel map.

Parameters
[in]box_min,box_maxThe minimum and maximum X/Y/Z unit coordinates that make up the box to be sampled
Returns
The average voxel map value of the sample box

Definition at line 979 of file voxel_map.cpp.

References _map_coord_idx_lower(), _map_coord_idx_upper(), _num_points, _sample_box_1d(), _sample_box_2d(), and _sample_box_3d().

◆ _sample_box_1d()

Real FEAT::Geometry::VoxelMap::_sample_box_1d ( i64  px0,
i64  px1,
i64  xidx0,
i64  xidx1,
i64  yidx = i64(0),
i64  zidx = i64(0) 
) const
protected

Helper function for _sample_box: samples a 1D edge in X line.

...the left-most interval equals to s0*(v0+V(idx0+1))/2 ...the right-most interval equals to s1*(v1+V(idx1-1))/2

Definition at line 750 of file voxel_map.cpp.

References _calc_sample_rate(), _check_point(), _num_points, _out_of_bounds_value, and ASSERT.

Referenced by _sample_box(), and _sample_box_2d().

◆ _sample_box_2d()

Real FEAT::Geometry::VoxelMap::_sample_box_2d ( i64  px0,
i64  px1,
i64  py0,
i64  py1,
i64  xidx0,
i64  xidx1,
i64  yidx0,
i64  yidx1,
i64  zidx = i64(0) 
) const
protected

Helper function for _sample_box: samples a 2D rectangle in XY plane.

Definition at line 843 of file voxel_map.cpp.

References _calc_sample_rate(), _num_points, _out_of_bounds_value, _sample_box_1d(), and ASSERT.

Referenced by _sample_box(), and _sample_box_3d().

◆ _sample_box_3d()

Real FEAT::Geometry::VoxelMap::_sample_box_3d ( i64  px0,
i64  px1,
i64  py0,
i64  py1,
i64  pz0,
i64  pz1,
i64  xidx0,
i64  xidx1,
i64  yidx0,
i64  yidx1,
i64  zidx0,
i64  zidx1 
) const
protected

Helper function for _sample_box: samples a 3D cuboid in XYZ volume.

Definition at line 900 of file voxel_map.cpp.

References _calc_sample_rate(), _num_points, _out_of_bounds_value, _sample_box_2d(), and ASSERT.

Referenced by _sample_box().

◆ _sample_point() [1/3]

Real FEAT::Geometry::VoxelMap::_sample_point ( const std::array< i64, 1 > &  p) const
protected

Samples the voxel map for a given 1D point by multi-linear interpolation.

Parameters
[in]pThe unit coordinates of the point whose voxel map value is to be sampled.
Returns
The multi-linearly interpolated value of the eight voxels surrounding the point

Definition at line 734 of file voxel_map.cpp.

References _map_coord_idx_lower(), and _sample_point_1d_x().

Referenced by sample_point().

◆ _sample_point() [2/3]

Real FEAT::Geometry::VoxelMap::_sample_point ( const std::array< i64, 2 > &  p) const
protected

Samples the voxel map for a given 2D point by multi-linear interpolation.

Parameters
[in]pThe unit coordinates of the point whose voxel map value is to be sampled.
Returns
The multi-linearly interpolated value of the eight voxels surrounding the point

Definition at line 739 of file voxel_map.cpp.

References _map_coord_idx_lower(), and _sample_point_2d().

◆ _sample_point() [3/3]

Real FEAT::Geometry::VoxelMap::_sample_point ( const std::array< i64, 3 > &  p) const
protected

Samples the voxel map for a given 3D point by multi-linear interpolation.

Parameters
[in]pThe unit coordinates of the point whose voxel map value is to be sampled.
Returns
The multi-linearly interpolated value of the eight voxels surrounding the point

Definition at line 744 of file voxel_map.cpp.

References _map_coord_idx_lower(), and _sample_point_3d().

◆ _sample_point_1d_x()

Real FEAT::Geometry::VoxelMap::_sample_point_1d_x ( i64  px,
i64  xidx,
i64  yidx = i64(0),
i64  zidx = i64(0) 
) const
protected

Helper function for _sample_point: samples a 1D edge in X-parallel line.

Definition at line 699 of file voxel_map.cpp.

References _calc_sample_rate(), and check_point().

Referenced by _sample_point(), and _sample_point_2d().

◆ _sample_point_1d_y()

Real FEAT::Geometry::VoxelMap::_sample_point_1d_y ( i64  py,
i64  xidx,
i64  yidx,
i64  zidx = i64(0) 
) const
protected

Helper function for _sample_point: samples a 1D edge in Y-parallel line.

Definition at line 706 of file voxel_map.cpp.

References _calc_sample_rate(), and check_point().

◆ _sample_point_1d_z()

Real FEAT::Geometry::VoxelMap::_sample_point_1d_z ( i64  pz,
i64  xidx,
i64  yidx,
i64  zidx 
) const
protected

Helper function for _sample_point: samples a 1D edge in Z-parallel line.

Definition at line 713 of file voxel_map.cpp.

References _calc_sample_rate(), and check_point().

◆ _sample_point_2d()

Real FEAT::Geometry::VoxelMap::_sample_point_2d ( i64  px,
i64  py,
i64  xidx,
i64  yidx,
i64  zidx = i64(0) 
) const
protected

Helper function for _sample_point: samples a 2D square in XY-parallel plane.

Definition at line 720 of file voxel_map.cpp.

References _calc_sample_rate(), and _sample_point_1d_x().

Referenced by _sample_point(), and _sample_point_3d().

◆ _sample_point_3d()

Real FEAT::Geometry::VoxelMap::_sample_point_3d ( i64  px,
i64  py,
i64  pz,
i64  xidx,
i64  yidx,
i64  zidx 
) const
protected

Helper function for _sample_point: samples a 3D cube in XYZ volume.

Definition at line 727 of file voxel_map.cpp.

References _calc_sample_rate(), and _sample_point_2d().

Referenced by _sample_point().

◆ calc_line_stride()

static u64 FEAT::Geometry::VoxelMap::calc_line_stride ( u64  num_x)
inlinestatic

ensure that we're not compiling on some freak architecture...

auxiliary function: compute default line stride from number of points in X dimension

Definition at line 448 of file voxel_map.hpp.

Referenced by set_num_points().

◆ check_box()

template<typename Coord_ , int dim_>
bool FEAT::Geometry::VoxelMap::check_box ( const Tiny::Vector< Coord_, dim_ > &  box_min,
const Tiny::Vector< Coord_, dim_ > &  box_max 
) const
inline

Checks the voxel map entries for a given box of voxels.

Parameters
[in]box_minThe minimum coordinates for the voxel box that is to be checked
[in]box_maxThe maximum coordinates for the voxel box that is to be checked
Note
If the box represented by box_min and box_max is clamped to the bounding box of the voxel map.
Returns
true, if at least one voxel contained in the box given by box_min and box_max is set to true, otherwise false.

Definition at line 1070 of file voxel_map.hpp.

References _check_box(), _map_coord_idx_lower(), _map_coord_idx_upper(), and unit_size.

◆ check_point()

bool FEAT::Geometry::VoxelMap::check_point ( u64  ix,
u64  iy = 0u,
u64  iz = 0u 
) const
inline

Returns the value of a specific voxel point in the voxel map.

Parameters
[in]ix,iy,izThe X/Y/Z indices of the point in the voxel map; each must be in range {0, ..., num_points[dim]-1}.
Returns
The value of the point at the desired voxel map position.

Definition at line 1006 of file voxel_map.hpp.

References _check_point(), _num_points, and ASSERT.

Referenced by _sample_point_1d_x(), _sample_point_1d_y(), and _sample_point_1d_z().

◆ check_point_nearest()

template<typename Coord_ , int dim_>
bool FEAT::Geometry::VoxelMap::check_point_nearest ( const Tiny::Vector< Coord_, dim_ > &  point) const
inline

Checks the voxel map entry for a given point.

Parameters
[in]pointThe point whose voxel map entry is to be checked.
Returns
The value of the voxel map in the position represented by point if the point is inside the voxel map bounding box, otherwise returns the out-of-bounds-value.

Definition at line 1025 of file voxel_map.hpp.

References _check_point_nearest(), and unit_size.

◆ compute_map()

template<typename Coord_ , int dim_>
void FEAT::Geometry::VoxelMap::compute_map ( const Dist::Comm comm,
VoxelMasker< Coord_, dim_ > &  masker,
bool  gather_to_all = true 
)
inline

Creates the voxel map based on a given masker object by utilizing MPI parallelization.

In an MPI-parallel simulation, this function must be called collectively by all processes in the communicator comm, since this function splits the workload across all these processes.

Attention
This function can only be called once the bounding box dimensions as well as the number of points or the resolution have been set.
Parameters
[in]commA transient reference to the communicator to use
[in]maskerA transient reference to the masker object
[in]gather_to_allSpecifies whether the voxel map is to be gathered on all processes in the communicator. The only reason to set this to false is when the map is only required on rank 0 for file output.

Definition at line 729 of file voxel_map.hpp.

References _compute_voxel_map(), _voxel_map, and XASSERT.

◆ compute_map_from_chart()

template<typename MeshType_ >
void FEAT::Geometry::VoxelMap::compute_map_from_chart ( const Dist::Comm comm,
const Geometry::Atlas::ChartBase< MeshType_ > &  chart,
bool  invert,
bool  gather_to_all = true 
)
inline

Creates the voxel map based on a chart by utilizing MPI parallelization.

In an MPI-parallel simulation, this function must be called collectively by all processes in the communicator comm, since this function splits the workload across all these processes.

Attention
This function can only be called once the bounding box dimensions as well as the number of points or the resolution have been set.
Parameters
[in]commA transient reference to the communicator to use
[in]chartA transient reference to the chart
[in]invertSpecifies whether the definition of inside and outside of the chart is to be swapped.
[in]gather_to_allSpecifies whether the voxel map is to be gathered on all processes in the communicator. The only reason to set this to false is when the map is only required on rank 0 for file output.

Definition at line 760 of file voxel_map.hpp.

References _compute_voxel_map(), _voxel_map, and XASSERT.

◆ compute_map_from_formula_2d()

void FEAT::Geometry::VoxelMap::compute_map_from_formula_2d ( const Dist::Comm comm,
const String formula,
bool  gather_to_all = true 
)

Creates the voxel map based on a 2D formula by utilizing MPI parallelization.

In an MPI-parallel simulation, this function must be called collectively by all processes in the communicator comm, since this function splits the workload across all these processes.

Attention
This function can only be called once the bounding box dimensions as well as the number of points or the resolution have been set.
Parameters
[in]commA transient reference to the communicator to use
[in]formulaA formula in x,y that evaluates to a value >= 0, if the points is inside the masked domain or to a value <= 0 otherwise.
[in]gather_to_allSpecifies whether the voxel map is to be gathered on all processes in the communicator. The only reason to set this to false is when the map is only required on rank 0 for file output.

Definition at line 136 of file voxel_map.cpp.

References _compute_voxel_map(), and XABORTM.

◆ compute_map_from_formula_3d()

void FEAT::Geometry::VoxelMap::compute_map_from_formula_3d ( const Dist::Comm comm,
const String formula,
bool  gather_to_all = true 
)

Creates the voxel map based on a 3D formula by utilizing MPI parallelization.

In an MPI-parallel simulation, this function must be called collectively by all processes in the communicator comm, since this function splits the workload across all these processes.

Attention
This function can only be called once the bounding box dimensions as well as the number of points or the resolution have been set.
Parameters
[in]commA transient reference to the communicator to use
[in]formulaA formula in x,y,z that evaluates to a value >= 0, if the points is inside the masked domain or to a value <= 0 otherwise.
[in]gather_to_allSpecifies whether the voxel map is to be gathered on all processes in the communicator. The only reason to set this to false is when the map is only required on rank 0 for file output.

Definition at line 149 of file voxel_map.cpp.

References _compute_voxel_map(), and XABORTM.

◆ compute_map_from_lambda_2d()

template<typename Lambda_ >
void FEAT::Geometry::VoxelMap::compute_map_from_lambda_2d ( const Dist::Comm comm,
Lambda_ &&  lambda,
bool  gather_to_all = true 
)
inline

Creates the voxel map based on a 2D lambda expression by utilizing MPI parallelization.

In an MPI-parallel simulation, this function must be called collectively by all processes in the communicator comm, since this function splits the workload across all these processes.

Attention
This function can only be called once the bounding box dimensions as well as the number of points or the resolution have been set.
Parameters
[in]commA transient reference to the communicator to use
[in]lambdaA lambda expression that accepts a Tiny::Vector<Real,2> as the only input parameter and returns true, if the point is inside the masked domain or false otherwise. This expression must be identical on all MPI processes.
[in]gather_to_allSpecifies whether the voxel map is to be gathered on all processes in the communicator. The only reason to set this to false is when the map is only required on rank 0 for file output.

Definition at line 791 of file voxel_map.hpp.

References _compute_voxel_map().

◆ compute_map_from_lambda_3d()

template<typename Lambda_ >
void FEAT::Geometry::VoxelMap::compute_map_from_lambda_3d ( const Dist::Comm comm,
Lambda_ &&  lambda,
bool  gather_to_all = true 
)
inline

Creates the voxel map based on a 3D lambda expression by utilizing MPI parallelization.

In an MPI-parallel simulation, this function must be called collectively by all processes in the communicator comm, since this function splits the workload across all these processes.

Attention
This function can only be called once the bounding box dimensions as well as the number of points or the resolution have been set.
Parameters
[in]commA transient reference to the communicator to use
[in]lambdaA lambda expression that accepts a Tiny::Vector<Real,3> as the only input parameter and returns true, if the point is inside the masked domain or false otherwise. This expression must be identical on all MPI processes.
[in]gather_to_allSpecifies whether the voxel map is to be gathered on all processes in the communicator. The only reason to set this to false is when the map is only required on rank 0 for file output.

Definition at line 821 of file voxel_map.hpp.

References _compute_voxel_map().

◆ compute_map_from_off_3d() [1/2]

void FEAT::Geometry::VoxelMap::compute_map_from_off_3d ( const Dist::Comm comm,
const String filename,
bool  invert,
bool  gather_to_all = true 
)

Creates the voxel map based on a 3D OFF model handled by CGAL by utilizing MPI parallelization.

In an MPI-parallel simulation, this function must be called collectively by all processes in the communicator comm, since this function splits the workload across all these processes.

Attention
This function can only be called once the bounding box dimensions as well as the number of points or the resolution have been set.
Parameters
[in]commA transient reference to the communicator to use
[in]filenameThe filename of the OFF model file that is to be read in. This file will be read in using the DistFileIO functionality and broadcasted across all processes in the communicator.
[in]invertSpecifies whether the definition of inside and outside of the model is to be swapped.
[in]gather_to_allSpecifies whether the voxel map is to be gathered on all processes in the communicator. The only reason to set this to false is when the map is only required on rank 0 for file output.

Definition at line 162 of file voxel_map.cpp.

References _bbox_min, compute_map_from_off_3d(), FEAT::DistFileIO::read_common(), and XASSERTM.

Referenced by compute_map_from_off_3d().

◆ compute_map_from_off_3d() [2/2]

void FEAT::Geometry::VoxelMap::compute_map_from_off_3d ( const Dist::Comm comm,
std::istream &  is,
bool  invert,
bool  gather_to_all = true 
)

Creates the voxel map based on a 3D OFF model handled by CGAL by utilizing MPI parallelization.

In an MPI-parallel simulation, this function must be called collectively by all processes in the communicator comm, since this function splits the workload across all these processes.

Attention
This function can only be called once the bounding box dimensions as well as the number of points or the resolution have been set.
Parameters
[in]commA transient reference to the communicator to use
[in]isThe input stream that the OFF file contents are to be parsed from. Must have the identical content on all MPI processes.
[in]Specifieswhether the definition of inside and outside of the model is to be swapped.
[in]gather_to_allSpecifies whether the voxel map is to be gathered on all processes in the communicator. The only reason to set this to false is when the map is only required on rank 0 for file output.

Definition at line 170 of file voxel_map.cpp.

References _bbox_min, _compute_voxel_map(), FEAT::Geometry::fm_off, XABORTM, and XASSERTM.

◆ export_plane_to_bmp() [1/2]

void FEAT::Geometry::VoxelMap::export_plane_to_bmp ( const String filename,
Real  z_coord 
) const

Exports a single plane of the voxel map to a monochrome BMP image file.

This function will export an entire plane of the voxel map in full resolution as a monochrome bitmap files.

Note
If you want to export the voxel map plane in a lower resolution, consider using the render_plane_to_bmp() function instead.
Parameters
[in]filenameThe filename for the BMP image file
[in]z_coordThe Z-coordinate of the plane that is to be exported.

Definition at line 469 of file voxel_map.cpp.

References export_plane_to_bmp().

Referenced by export_plane_to_bmp().

◆ export_plane_to_bmp() [2/2]

void FEAT::Geometry::VoxelMap::export_plane_to_bmp ( std::ostream &  os,
Real  z_coord 
) const

Exports a single plane of the voxel map to a monochrome BMP image file.

This function will export an entire plane of the voxel map in full resolution as a monochrome bitmap files.

Note
If you want to export the voxel map plane in a lower resolution, consider using the render_plane_to_bmp() function instead.
Parameters
[out]osThe binary output stream to write to
[in]z_coordThe Z-coordinate of the plane that is to be exported.

Definition at line 476 of file voxel_map.cpp.

References _export_plane_to_bmp(), _map_coord_idx_nearest(), and unit_size.

◆ export_to_bmp()

void FEAT::Geometry::VoxelMap::export_to_bmp ( const String filename_prefix) const

Exports the voxel map to a sequence of monochrome BMP image files.

This function will export the entire voxel map in full resolution as a sequence of monochrome bitmap files. The maximum voxel map dimension is 20Kx20Kx20K, which would result in 20000 BMP files with a resolution 400 MPixels each.

Note
If you want to export the voxel map in a lower resolution, consider using the render_to_bmp() function instead.
Parameters
[in]filename_prefixThe filename prefix for the image file sequence.

Definition at line 456 of file voxel_map.cpp.

References _export_plane_to_bmp(), _num_planes, _num_points, FEAT::String::pad_front(), FEAT::stringify(), and XASSERTM.

◆ get_bounding_box_max()

Real FEAT::Geometry::VoxelMap::get_bounding_box_max ( int  dim) const
inline

Returns the maximum bounding box coordinate for a given dimensions.

Definition at line 601 of file voxel_map.hpp.

References unit_size, and XASSERT.

◆ get_bounding_box_min()

Real FEAT::Geometry::VoxelMap::get_bounding_box_min ( int  dim) const
inline

Returns the minimum bounding box coordinate for a given dimensions.

Definition at line 594 of file voxel_map.hpp.

References unit_size, and XASSERT.

◆ get_bounding_box_volume()

Real FEAT::Geometry::VoxelMap::get_bounding_box_volume ( ) const

Returns the volume of the 1D/2D/3D voxel map bounding box.

Definition at line 507 of file voxel_map.cpp.

References _bbox_min, _num_points, FEAT::Math::cub(), FEAT::Math::sqr(), and unit_size.

Referenced by get_domain_volume().

◆ get_domain_coverage()

Real FEAT::Geometry::VoxelMap::get_domain_coverage ( ) const
inline

Returns the domain coverage of the voxel map.

The domain coverage of the voxel map is defined as the number of voxels of value 1 divided by the total number of voxels. The returned value is always in range [0,1].

Returns
The domain coverage of the voxel map

Definition at line 695 of file voxel_map.hpp.

References unit_size.

Referenced by get_domain_volume().

◆ get_domain_volume()

Real FEAT::Geometry::VoxelMap::get_domain_volume ( ) const
inline

Returns the covered domain volume.

Definition at line 702 of file voxel_map.hpp.

References get_bounding_box_volume(), and get_domain_coverage().

◆ get_map()

const std::vector< char > & FEAT::Geometry::VoxelMap::get_map ( ) const
inline

Returns a reference to the voxel map.

Definition at line 682 of file voxel_map.hpp.

References _voxel_map.

◆ get_num_points()

Index FEAT::Geometry::VoxelMap::get_num_points ( int  dim) const
inline

Returns the number of points of the voxel map in each dimension.

Definition at line 641 of file voxel_map.hpp.

◆ get_num_voxels()

Index FEAT::Geometry::VoxelMap::get_num_voxels ( ) const
inline

Returns the total number of voxels in the voxel map.

Definition at line 658 of file voxel_map.hpp.

Referenced by _compute_domain_coverage().

◆ get_out_of_bounds_value()

bool FEAT::Geometry::VoxelMap::get_out_of_bounds_value ( ) const
inline

Returns the out-of-bounds value for the voxel map.

Definition at line 627 of file voxel_map.hpp.

References _out_of_bounds_value.

◆ get_stride_line()

u64 FEAT::Geometry::VoxelMap::get_stride_line ( ) const
inline

Returns the line stride, i.e. the size of a single X-line in bytes.

Definition at line 664 of file voxel_map.hpp.

References _stride_line.

◆ get_stride_plane()

u64 FEAT::Geometry::VoxelMap::get_stride_plane ( ) const
inline

Returns the plane stride, i.e. the size of a single XY-plane in bytes.

Definition at line 670 of file voxel_map.hpp.

References _stride_plane.

◆ get_stride_volume()

u64 FEAT::Geometry::VoxelMap::get_stride_volume ( ) const
inline

Returns the volume stride, i.e. the size of a single XYZ-volume in bytes.

Definition at line 676 of file voxel_map.hpp.

◆ operator=()

VoxelMap & FEAT::Geometry::VoxelMap::operator= ( VoxelMap &&  other)

◆ read() [1/2]

VoxelMap::ReadResult FEAT::Geometry::VoxelMap::read ( const Dist::Comm comm,
const String filename 
)

Reads a voxel map from a file.

Parameters
[in]commA transient reference to a communicator.
[in]filenameThe filename of the voxel map file to be read in.
Returns
A ReadResult object that contains more detailed information about the input file, including bytes read and total compression buffer size. Can be cast to bool to check whether the read was successful.

Definition at line 322 of file voxel_map.cpp.

References read(), and FEAT::DistFileIO::read_common().

Referenced by read().

◆ read() [2/2]

VoxelMap::ReadResult FEAT::Geometry::VoxelMap::read ( std::istream &  is,
const String filename = "" 
)

Reads a voxel map from a binary input stream.

Parameters
[in]isA transient reference to a binary input stream to read form
[in]filenameThe filename of the voxel map file to be read in. This is just used for error output.

Definition at line 329 of file voxel_map.cpp.

References _bbox_min, _coverage, _num_lines, _num_planes, _num_points, _stride_line, _stride_plane, _voxel_map, header_size, magic_number, FEAT::Math::min(), FEAT::stringify(), XABORTM, and XASSERTM.

◆ render_plane_to_bmp() [1/2]

void FEAT::Geometry::VoxelMap::render_plane_to_bmp ( const String filename,
Index  width,
Index  height,
Real  z_min,
Real  z_max 
) const

Renders a plane range of voxel map to a single gray-scale BMP image file.

This function will render a range of plane of the voxel map in reduced resolution as a single gray-scale bitmap file, where each pixel contains the average value of all voxels that have been reduced to it.

Parameters
[in]filenameThe filename for the BMP image
[in]widthThe desired bitmap width
[in]heightThe desired bitmap height
[in]z_min,z_maxThe Z-coordinates of the plane range that has to be sampled into a single output plane

Definition at line 495 of file voxel_map.cpp.

References render_plane_to_bmp().

Referenced by render_plane_to_bmp().

◆ render_plane_to_bmp() [2/2]

void FEAT::Geometry::VoxelMap::render_plane_to_bmp ( std::ostream &  os,
Index  width,
Index  height,
Real  z_min,
Real  z_max 
) const

Renders a plane range of voxel map to a single gray-scale BMP image file stream.

This function will render a range of plane of the voxel map in reduced resolution as a single gray-scale bitmap file, where each pixel contains the average value of all voxels that have been reduced to it.

Parameters
[out]osThe output stream to write to
[in]widthThe desired bitmap width
[in]heightThe desired bitmap height
[in]z_min,z_maxThe Z-coordinates of the plane range that has to be sampled into a single output plane

Definition at line 502 of file voxel_map.cpp.

References _render_plane_to_bmp(), and unit_size.

◆ render_to_bmp()

void FEAT::Geometry::VoxelMap::render_to_bmp ( const String filename_prefix,
Index  width,
Index  height,
Index  depth 
) const

Renders the voxel map to a sequence of gray-scale BMP image files.

This function will render the entire voxel map in reduced resolution as a sequence of gray-scale bitmap files, where each pixel contains the average value of all voxels that have been reduced to it.

Parameters
[in]filename_prefixThe filename prefix for the image file sequence
[in]widthThe desired bitmap width
[in]heightThe desired bitmap height
[in]depthThe desired number of bitmap images in the sequence

Definition at line 481 of file voxel_map.cpp.

References _bbox_min, _render_plane_to_bmp(), FEAT::String::pad_front(), FEAT::stringify(), and XASSERTM.

◆ sample_box()

template<typename Coord_ , int dim_>
Real FEAT::Geometry::VoxelMap::sample_box ( const Tiny::Vector< Coord_, dim_ > &  box_min,
const Tiny::Vector< Coord_, dim_ > &  box_max 
) const
inline

Samples the voxel map entries for a given box of voxels.

Parameters
[in]box_minThe minimum coordinates for the voxel box that is to be sampled
[in]box_maxThe maximum coordinates for the voxel box that is to be sampled
Note
The box represented by box_min and box_max is clamped to the bounding box of the voxel map.
Returns
The average number of voxel map entries in the box which are set to true.
Note
The return value is 1 if all voxels in the box are set to 1 and the return value is 0 if all voxels in the box are 0; in any other case the return value is > 0 and < 1, representing the average value of all voxels in the sampled box.

Definition at line 1122 of file voxel_map.hpp.

References _sample_box(), and unit_size.

◆ sample_point()

template<typename Coord_ , int dim_>
Real FEAT::Geometry::VoxelMap::sample_point ( const Tiny::Vector< Coord_, dim_ > &  point) const
inline

Samples the voxel map entry for a given point by multi-linear interpolation.

Parameters
[in]pointThe point whose voxel map entry is to be samples.
Returns
The multi-linearly interpolated value of the voxel map in the position represented by point.

Definition at line 1044 of file voxel_map.hpp.

References _sample_point(), and unit_size.

◆ set_bounding_box_2d()

void FEAT::Geometry::VoxelMap::set_bounding_box_2d ( Real  x_min,
Real  x_max,
Real  y_min,
Real  y_max 
)

Sets the bounding box for a 2D voxel map.

Parameters
[in]x_min,x_maxThe range of the voxel map in X dimension
[in]y_min,y_maxThe range of the voxel map in Y dimension

Definition at line 79 of file voxel_map.cpp.

References _bbox_min, unit_size, and XASSERTM.

◆ set_bounding_box_3d()

void FEAT::Geometry::VoxelMap::set_bounding_box_3d ( Real  x_min,
Real  x_max,
Real  y_min,
Real  y_max,
Real  z_min,
Real  z_max 
)

Sets the bounding box for a 3D voxel map.

Parameters
[in]x_min,x_maxThe range of the voxel map in X dimension
[in]y_min,y_maxThe range of the voxel map in Y dimension
[in]z_min,z_maxThe range of the voxel map in Z dimension

Definition at line 92 of file voxel_map.cpp.

References _bbox_min, unit_size, and XASSERTM.

◆ set_num_points()

void FEAT::Geometry::VoxelMap::set_num_points ( Index  num_x,
Index  num_y,
Index  num_z = 0u 
)

Sets the number of points of the voxel map in each dimension.

Parameters
[in]num_x,num_y,num_zThe number of points in each dimension

Definition at line 108 of file voxel_map.cpp.

References _bbox_min, _num_lines, _num_planes, _num_points, _stride_line, _stride_plane, _voxel_map, calc_line_stride(), and XASSERT.

Referenced by set_resolution().

◆ set_out_of_bounds_value()

void FEAT::Geometry::VoxelMap::set_out_of_bounds_value ( bool  value)
inline

Sets the out-of-bounds value for the voxel map.

The out-of-bounds-value specifies whether the outside of the voxel map bounding box is being treated as 'inside' (true) or 'outside' (false), which is the value that is returned by the check_point and check_box functions if the point/bounding box is not inside of the voxel map bounding box.

Parameters
[in]valueThe out-of-bounds-value for the voxel map

Definition at line 621 of file voxel_map.hpp.

References FEAT::value.

◆ set_resolution()

void FEAT::Geometry::VoxelMap::set_resolution ( Real  max_res)

Sets the resolution for the voxel map, i.e. the maximum voxel size in each dimension.

This function computes the number of points in each dimension to ensure that the size of a voxel is less or equal to max_res in any dimension.

Parameters
[in]max_resThe maximum voxel size; must be > 0.

Definition at line 127 of file voxel_map.cpp.

References _bbox_min, set_num_points(), unit_size, and XASSERTM.

◆ write() [1/2]

VoxelMap::WriteResult FEAT::Geometry::VoxelMap::write ( const String filename,
const u64  compress_block_size = 128u,
const int  compress_level = 9 
) const

Writes the voxel map into a file that can later be read in.

Attention
In an MPI-parallel simulation, only one process shall execute this function!
Parameters
[in]filenameThe filename of the output file. The recommended extension is 'vxl'.
[in]compress_block_sizeThe maximum size of a single ZLIB compression block in MB. Set to 0 to disable ZLIB compression.
[in]compress_levelThe ZLIB compression level in range [1,9], where 1 is fastest compression and 9 is best compression.
Returns
A WriteResult object that contains more detailed information about the output file, including bytes written and total compression buffer size. Can be cast to bool to check whether the write was successful.

Definition at line 186 of file voxel_map.cpp.

References write().

Referenced by write().

◆ write() [2/2]

VoxelMap::WriteResult FEAT::Geometry::VoxelMap::write ( std::ostream &  os,
const u64  compress_block_size = 128u,
const int  compress_level = 9 
) const

Writes the voxel map into a binary output stream.

Attention
In an MPI-parallel simulation, only one process shall execute this function!
Parameters
[out]osThe binary output stream to write to.
[in]compress_block_sizeThe maximum size of a single ZLIB compression block in MB. Set to 0 to disable ZLIB compression.
[in]compress_levelThe ZLIB compression level in range [1,9], where 1 is fastest compression and 9 is best compression.
Returns
A WriteResult object that contains more detailed information about the output file, including bytes written and total compression buffer size. Can be cast to bool to check whether the write was successful.

Definition at line 195 of file voxel_map.cpp.

References _bbox_min, _coverage, _num_planes, _num_points, _stride_line, _stride_plane, _voxel_map, magic_number, FEAT::Math::max(), FEAT::Math::min(), XABORTM, XASSERT, and XASSERTM.

Member Data Documentation

◆ _bbox_max

std::array<i64, 3u> FEAT::Geometry::VoxelMap::_bbox_max
protected

Definition at line 533 of file voxel_map.hpp.

◆ _bbox_min

◆ _coverage

u64 FEAT::Geometry::VoxelMap::_coverage
protected

the coverage of the domain, relative to unit size

Definition at line 549 of file voxel_map.hpp.

Referenced by operator=(), read(), and write().

◆ _create_stage

int FEAT::Geometry::VoxelMap::_create_stage
protected

current creation stage

Definition at line 531 of file voxel_map.hpp.

Referenced by operator=().

◆ _num_lines

u64 FEAT::Geometry::VoxelMap::_num_lines
protected

number of lines; 2D = num_points[2]; 3D: = num_points[2] * num_points[3]

Definition at line 541 of file voxel_map.hpp.

Referenced by _compress_voxel_map_line(), operator=(), read(), and set_num_points().

◆ _num_planes

u64 FEAT::Geometry::VoxelMap::_num_planes
protected

number of planes; 2D = 1; 3D = num_points[3]

Definition at line 543 of file voxel_map.hpp.

Referenced by export_to_bmp(), operator=(), read(), set_num_points(), and write().

◆ _num_points

◆ _out_of_bounds_value

bool FEAT::Geometry::VoxelMap::_out_of_bounds_value
protected

the out-of-bounds-value for the voxel map

Definition at line 547 of file voxel_map.hpp.

Referenced by _check_point(), _check_point_nearest(), _sample_box_1d(), _sample_box_2d(), _sample_box_3d(), get_out_of_bounds_value(), and operator=().

◆ _stride_line

u64 FEAT::Geometry::VoxelMap::_stride_line
protected

stride of a single voxel line in X dimension

Definition at line 537 of file voxel_map.hpp.

Referenced by _check_box(), _check_point(), _check_point_nearest(), _compress_voxel_map_line(), _export_plane_to_bmp(), get_stride_line(), operator=(), read(), set_num_points(), and write().

◆ _stride_plane

u64 FEAT::Geometry::VoxelMap::_stride_plane
protected

stride of a single plane in XY dimensions

Definition at line 539 of file voxel_map.hpp.

Referenced by _check_point(), _export_plane_to_bmp(), get_stride_plane(), operator=(), read(), set_num_points(), and write().

◆ _voxel_map

std::vector<char> FEAT::Geometry::VoxelMap::_voxel_map
protected

the actual voxel map; its size may be larger than necessary due to padding

Definition at line 545 of file voxel_map.hpp.

Referenced by _check_box(), _check_point(), _check_point_nearest(), _compress_voxel_map_line(), _compute_domain_coverage(), _export_plane_to_bmp(), compute_map(), compute_map_from_chart(), get_map(), operator=(), read(), set_num_points(), and write().

◆ header_size

constexpr u64 FEAT::Geometry::VoxelMap::header_size = 136ull
staticconstexpr

size of voxel map header in bytes

Definition at line 417 of file voxel_map.hpp.

Referenced by read().

◆ magic_number

constexpr u64 FEAT::Geometry::VoxelMap::magic_number = 0x70614D6C65786F56ull
staticconstexpr

voxel map magic number

Definition at line 414 of file voxel_map.hpp.

Referenced by read(), and write().

◆ unit_size

constexpr i64 FEAT::Geometry::VoxelMap::unit_size = 1'000'000'000
staticconstexpr

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