SemiDiscreteOT 1.0
Semi-Discrete Optimal Transport Library
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Attributes | List of all members
VTKHandler< dim, spacedim > Class Template Reference

A handler class for reading and interpolating data from VTK files. More...

#include <VtkHandler.h>

Inheritance diagram for VTKHandler< dim, spacedim >:
Inheritance graph
[legend]
Collaboration diagram for VTKHandler< dim, spacedim >:
Collaboration graph
[legend]

Public Types

enum class  DataLocation { PointData , CellData }
 Enum to specify whether the data is associated with points or cells. More...
 

Public Member Functions

 VTKHandler (const std::string &filename, const bool is_binary=false, const double scaling_factor=1.0)
 Constructor for the VTKHandler.
 
void read_file ()
 Reads the VTK file and initializes the internal data structures.
 
void setup_field (const std::string &field_name, const DataLocation data_location=DataLocation::PointData, const unsigned int component=0)
 Sets up the field to be used for interpolation.
 
virtual double value (const Point< spacedim > &p, const unsigned int component=0) const override
 Interpolates the value of the selected field at a given point.
 
vtkSmartPointer< vtkUnstructuredGrid > get_grid () const
 Returns the underlying VTK unstructured grid.
 
vtkSmartPointer< vtkDataArray > get_field_data () const
 Returns the data array for the selected field.
 
int get_num_components () const
 Returns the number of components in the selected field data.
 

Private Attributes

std::string filename
 
bool is_binary
 
double scaling_factor
 
vtkSmartPointer< vtkUnstructuredGrid > vtk_grid
 
vtkSmartPointer< vtkDataArray > field_data
 
vtkSmartPointer< vtkPointLocator > point_locator
 
vtkSmartPointer< vtkCellLocator > cell_locator
 
DataLocation data_location
 
std::string field_name
 
unsigned int selected_component
 

Detailed Description

template<int dim, int spacedim = dim>
class VTKHandler< dim, spacedim >

A handler class for reading and interpolating data from VTK files.

This class inherits from dealii::Function and provides functionality to read VTK unstructured grids, extract scalar or vector fields from point or cell data, and interpolate these fields at arbitrary points in space.

Template Parameters
dimThe dimension of the mesh.
spacedimThe dimension of the space the mesh is embedded in.

Definition at line 33 of file VtkHandler.h.

Member Enumeration Documentation

◆ DataLocation

template<int dim, int spacedim = dim>
enum class VTKHandler::DataLocation
strong

Enum to specify whether the data is associated with points or cells.

Enumerator
PointData 

Data is associated with the vertices of the grid.

CellData 

Data is associated with the cells of the grid.

Definition at line 39 of file VtkHandler.h.

Constructor & Destructor Documentation

◆ VTKHandler()

template<int dim, int spacedim>
VTKHandler< dim, spacedim >::VTKHandler ( const std::string &  filename,
const bool  is_binary = false,
const double  scaling_factor = 1.0 
)

Constructor for the VTKHandler.

Parameters
filenameThe path to the VTK file.
is_binaryWhether the VTK file is in binary format.
scaling_factorA factor to scale the data by upon reading.

Definition at line 15 of file VtkHandler.cc.

Here is the call graph for this function:

Member Function Documentation

◆ read_file()

template<int dim, int spacedim>
void VTKHandler< dim, spacedim >::read_file ( )

Reads the VTK file and initializes the internal data structures.

Exceptions
std::runtime_errorif the file cannot be opened or is invalid.

Definition at line 33 of file VtkHandler.cc.

Here is the caller graph for this function:

◆ setup_field()

template<int dim, int spacedim>
void VTKHandler< dim, spacedim >::setup_field ( const std::string &  field_name,
const DataLocation  data_location = DataLocation::PointData,
const unsigned int  component = 0 
)

Sets up the field to be used for interpolation.

Parameters
field_nameThe name of the data array in the VTK file.
data_locationWhether the data is point or cell data.
componentThe component of the data array to use (for vector fields).

Definition at line 94 of file VtkHandler.cc.

Here is the caller graph for this function:

◆ value()

template<int dim, int spacedim>
double VTKHandler< dim, spacedim >::value ( const Point< spacedim > &  p,
const unsigned int  component = 0 
) const
overridevirtual

Interpolates the value of the selected field at a given point.

This method overrides the value method of the dealii::Function base class.

Parameters
pThe point at which to interpolate the value.
componentThe component of the function to evaluate (unused, as the component is pre-selected in setup_field).
Returns
The interpolated value of the field at point p.

Definition at line 134 of file VtkHandler.cc.

◆ get_grid()

template<int dim, int spacedim = dim>
vtkSmartPointer< vtkUnstructuredGrid > VTKHandler< dim, spacedim >::get_grid ( ) const
inline

Returns the underlying VTK unstructured grid.

Returns
A vtkSmartPointer to the vtkUnstructuredGrid.

Definition at line 88 of file VtkHandler.h.

◆ get_field_data()

template<int dim, int spacedim = dim>
vtkSmartPointer< vtkDataArray > VTKHandler< dim, spacedim >::get_field_data ( ) const
inline

Returns the data array for the selected field.

Returns
A vtkSmartPointer to the vtkDataArray.

Definition at line 94 of file VtkHandler.h.

◆ get_num_components()

template<int dim, int spacedim = dim>
int VTKHandler< dim, spacedim >::get_num_components ( ) const
inline

Returns the number of components in the selected field data.

Returns
The number of components.

Definition at line 100 of file VtkHandler.h.

Member Data Documentation

◆ filename

template<int dim, int spacedim = dim>
std::string VTKHandler< dim, spacedim >::filename
private

Definition at line 103 of file VtkHandler.h.

◆ is_binary

template<int dim, int spacedim = dim>
bool VTKHandler< dim, spacedim >::is_binary
private

Definition at line 104 of file VtkHandler.h.

◆ scaling_factor

template<int dim, int spacedim = dim>
double VTKHandler< dim, spacedim >::scaling_factor
private

Definition at line 105 of file VtkHandler.h.

◆ vtk_grid

template<int dim, int spacedim = dim>
vtkSmartPointer<vtkUnstructuredGrid> VTKHandler< dim, spacedim >::vtk_grid
private

Definition at line 107 of file VtkHandler.h.

◆ field_data

template<int dim, int spacedim = dim>
vtkSmartPointer<vtkDataArray> VTKHandler< dim, spacedim >::field_data
private

Definition at line 108 of file VtkHandler.h.

◆ point_locator

template<int dim, int spacedim = dim>
vtkSmartPointer<vtkPointLocator> VTKHandler< dim, spacedim >::point_locator
private

Definition at line 109 of file VtkHandler.h.

◆ cell_locator

template<int dim, int spacedim = dim>
vtkSmartPointer<vtkCellLocator> VTKHandler< dim, spacedim >::cell_locator
private

Definition at line 110 of file VtkHandler.h.

◆ data_location

template<int dim, int spacedim = dim>
DataLocation VTKHandler< dim, spacedim >::data_location
private

Definition at line 112 of file VtkHandler.h.

◆ field_name

template<int dim, int spacedim = dim>
std::string VTKHandler< dim, spacedim >::field_name
private

Definition at line 113 of file VtkHandler.h.

◆ selected_component

template<int dim, int spacedim = dim>
unsigned int VTKHandler< dim, spacedim >::selected_component
private

Definition at line 114 of file VtkHandler.h.


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