Reference documentation for deal.II version 9.5.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Member Functions | Static Private Member Functions | Private Attributes | List of all members
GridTools::MarchingCubeAlgorithm< dim, VectorType > Class Template Reference

#include <deal.II/grid/grid_tools.h>

Public Types

using value_type = typename VectorType::value_type
 

Public Member Functions

 MarchingCubeAlgorithm (const Mapping< dim, dim > &mapping, const FiniteElement< dim, dim > &fe, const unsigned int n_subdivisions=1, const double tolerance=1e-10)
 
void process (const DoFHandler< dim > &background_dof_handler, const VectorType &ls_vector, const double iso_level, std::vector< Point< dim > > &vertices, std::vector< CellData< dim==1 ? 1 :dim - 1 > > &cells) const
 
void process (const DoFHandler< dim > &background_dof_handler, const VectorType &ls_vector, const double iso_level, std::vector< Point< dim > > &vertices) const
 
void process_cell (const typename DoFHandler< dim >::active_cell_iterator &cell, const VectorType &ls_vector, const double iso_level, std::vector< Point< dim > > &vertices, std::vector< CellData< dim==1 ? 1 :dim - 1 > > &cells) const
 
void process_cell (const typename DoFHandler< dim >::active_cell_iterator &cell, const VectorType &ls_vector, const double iso_level, std::vector< Point< dim > > &vertices) const
 

Private Member Functions

void process_cell (std::vector< value_type > &ls_values, const std::vector< Point< dim > > &points, const double iso_level, std::vector< Point< dim > > &vertices, std::vector< CellData< dim==1 ? 1 :dim - 1 > > &cells, const bool write_back_cell_data=true) const
 
void process_sub_cell (const std::vector< value_type > &, const std::vector< Point< 1 > > &, const std::vector< unsigned int > &, const double, std::vector< Point< 1 > > &, std::vector< CellData< 1 > > &, const bool) const
 
void process_sub_cell (const std::vector< value_type > &ls_values, const std::vector< Point< 2 > > &points, const std::vector< unsigned int > &mask, const double iso_level, std::vector< Point< 2 > > &vertices, std::vector< CellData< 1 > > &cells, const bool write_back_cell_data) const
 
void process_sub_cell (const std::vector< value_type > &ls_values, const std::vector< Point< 3 > > &points, const std::vector< unsigned int > &mask, const double iso_level, std::vector< Point< 3 > > &vertices, std::vector< CellData< 2 > > &cells, const bool write_back_cell_data) const
 

Static Private Member Functions

static Quadrature< dim > create_quadrature_rule (const unsigned int n_subdivisions)
 

Private Attributes

const unsigned int n_subdivisions
 
const double tolerance
 
FEValues< dim > fe_values
 

Detailed Description

template<int dim, typename VectorType>
class GridTools::MarchingCubeAlgorithm< dim, VectorType >

An implementation of the marching-square (2d) and marching-cube algorithm for creating data structures (vectors of Point and CellData) to create a linear/bilinear surface mesh on the iso line/contour of a scalar field.

To improve the approximation of the iso line/contour and the resulting linear surface mesh, one increases the number of subdivision so that the algorithm is not run on a cell but on subcells with vertex values having been interpolated from the cell values.

Note
The resulting mesh will contain lines in 2d and triangles in 3d.
The resulting mesh will not be of high quality, since it might contain cells with very small diameters if the mesh is cut close to a vertex.
Iso lines/contours as a saddle point within a subcell is not detected by the implemented algorithm.

Definition at line 3739 of file grid_tools.h.

Member Typedef Documentation

◆ value_type

template<int dim, typename VectorType >
using GridTools::MarchingCubeAlgorithm< dim, VectorType >::value_type = typename VectorType::value_type

Value type of vector.

Definition at line 3745 of file grid_tools.h.

Constructor & Destructor Documentation

◆ MarchingCubeAlgorithm()

template<int dim, typename VectorType >
GridTools::MarchingCubeAlgorithm< dim, VectorType >::MarchingCubeAlgorithm ( const Mapping< dim, dim > &  mapping,
const FiniteElement< dim, dim > &  fe,
const unsigned int  n_subdivisions = 1,
const double  tolerance = 1e-10 
)

Constructor.

Definition at line 7315 of file grid_tools.cc.

Member Function Documentation

◆ process() [1/2]

template<int dim, typename VectorType >
void GridTools::MarchingCubeAlgorithm< dim, VectorType >::process ( const DoFHandler< dim > &  background_dof_handler,
const VectorType &  ls_vector,
const double  iso_level,
std::vector< Point< dim > > &  vertices,
std::vector< CellData< dim==1 ? 1 :dim - 1 > > &  cells 
) const

Process all locally-owned cells and fill vertices and cells for all cells that are cut.

Note
This function is only implemented for dim>1. Use process(background_dof_handler, ls_vector, iso_level, vertices) for dim==1.
Duplicate vertices are not deleted.

Definition at line 7367 of file grid_tools.cc.

◆ process() [2/2]

template<int dim, typename VectorType >
void GridTools::MarchingCubeAlgorithm< dim, VectorType >::process ( const DoFHandler< dim > &  background_dof_handler,
const VectorType &  ls_vector,
const double  iso_level,
std::vector< Point< dim > > &  vertices 
) const

Process all locally-owned cells and fill vertices for all cells that are cut.

Definition at line 7387 of file grid_tools.cc.

◆ process_cell() [1/3]

template<int dim, typename VectorType >
void GridTools::MarchingCubeAlgorithm< dim, VectorType >::process_cell ( const typename DoFHandler< dim >::active_cell_iterator &  cell,
const VectorType &  ls_vector,
const double  iso_level,
std::vector< Point< dim > > &  vertices,
std::vector< CellData< dim==1 ? 1 :dim - 1 > > &  cells 
) const

Process the provided cell and fill vertices and cells for all cells that are cut.

Note
The resulting vectors are empty if the cell is not cut.
This function is only implemented for dim>1. Use process_cell(cell, ls_vector, iso_level, vertices) for dim==1.

Definition at line 7403 of file grid_tools.cc.

◆ process_cell() [2/3]

template<int dim, typename VectorType >
void GridTools::MarchingCubeAlgorithm< dim, VectorType >::process_cell ( const typename DoFHandler< dim >::active_cell_iterator &  cell,
const VectorType &  ls_vector,
const double  iso_level,
std::vector< Point< dim > > &  vertices 
) const

Process the provided cell and fill vertices for all cells that are cut.

Note
The resulting vector is empty if the cell is not cut.

Definition at line 7427 of file grid_tools.cc.

◆ create_quadrature_rule()

template<int dim, typename VectorType >
Quadrature< dim > GridTools::MarchingCubeAlgorithm< dim, VectorType >::create_quadrature_rule ( const unsigned int  n_subdivisions)
staticprivate

Internal function to create a quadrature rule with n_subdivisions+1 equally-positioned quadrature points.

Definition at line 7332 of file grid_tools.cc.

◆ process_cell() [3/3]

template<int dim, typename VectorType >
void GridTools::MarchingCubeAlgorithm< dim, VectorType >::process_cell ( std::vector< value_type > &  ls_values,
const std::vector< Point< dim > > &  points,
const double  iso_level,
std::vector< Point< dim > > &  vertices,
std::vector< CellData< dim==1 ? 1 :dim - 1 > > &  cells,
const bool  write_back_cell_data = true 
) const
private

Process a cell.

Definition at line 7453 of file grid_tools.cc.

◆ process_sub_cell() [1/3]

template<int dim, typename VectorType >
void GridTools::MarchingCubeAlgorithm< dim, VectorType >::process_sub_cell ( const std::vector< value_type > &  ,
const std::vector< Point< 1 > > &  ,
const std::vector< unsigned int > &  ,
const double  ,
std::vector< Point< 1 > > &  ,
std::vector< CellData< 1 > > &  ,
const bool   
) const
inlineprivate

Dummy function for 1d processing a sub-cell.

Definition at line 3831 of file grid_tools.h.

◆ process_sub_cell() [2/3]

template<int dim, typename VectorType >
void GridTools::MarchingCubeAlgorithm< dim, VectorType >::process_sub_cell ( const std::vector< value_type > &  ls_values,
const std::vector< Point< 2 > > &  points,
const std::vector< unsigned int > &  mask,
const double  iso_level,
std::vector< Point< 2 > > &  vertices,
std::vector< CellData< 1 > > &  cells,
const bool  write_back_cell_data 
) const
private

Process a sub-cell (2d).

Note
Subcells with saddle points are ignored. Please increase the number of subdivisions in this case.

Definition at line 7543 of file grid_tools.cc.

◆ process_sub_cell() [3/3]

template<int dim, typename VectorType >
void GridTools::MarchingCubeAlgorithm< dim, VectorType >::process_sub_cell ( const std::vector< value_type > &  ls_values,
const std::vector< Point< 3 > > &  points,
const std::vector< unsigned int > &  mask,
const double  iso_level,
std::vector< Point< 3 > > &  vertices,
std::vector< CellData< 2 > > &  cells,
const bool  write_back_cell_data 
) const
private

Process a sub-cell (3d).

Definition at line 7626 of file grid_tools.cc.

Member Data Documentation

◆ n_subdivisions

template<int dim, typename VectorType >
const unsigned int GridTools::MarchingCubeAlgorithm< dim, VectorType >::n_subdivisions
private

Number of subdivisions each cell is subdivided into in each direction to improve the approximation.

Definition at line 3873 of file grid_tools.h.

◆ tolerance

template<int dim, typename VectorType >
const double GridTools::MarchingCubeAlgorithm< dim, VectorType >::tolerance
private

Absolute tolerance specifying the minimum distance between a vertex and the cut point so that a line is considered cut.

Definition at line 3879 of file grid_tools.h.

◆ fe_values

template<int dim, typename VectorType >
FEValues<dim> GridTools::MarchingCubeAlgorithm< dim, VectorType >::fe_values
mutableprivate

FEValues used internally and set up with a quadrature rule with the correct number of subdivisions.

Definition at line 3885 of file grid_tools.h.


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