VTK  9.4.20241117
Functions
vtkConduitToDataObject Namespace Reference

Functions

VTKIOCATALYSTCONDUIT_EXPORT bool FillPartionedDataSet (vtkPartitionedDataSet *output, const conduit_cpp::Node &meshNode)
 Fill the vtkPartitionedDataSet input.
 
VTKIOCATALYSTCONDUIT_EXPORT bool FillAMRMesh (vtkOverlappingAMR *amr, const conduit_cpp::Node &node)
 Fill the vtkOverlappingAMR input.
 
VTKIOCATALYSTCONDUIT_EXPORT bool AddFieldData (vtkDataObject *output, const conduit_cpp::Node &stateFields, bool isAMReX=false)
 Add FieldData arrays to output data object.
 
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkPointsCreatePoints (const conduit_cpp::Node &coords)
 Create a vtkPoints from a coordset node that respect the following requirements:
 
VTKIOCATALYSTCONDUIT_EXPORT void SetPolyhedralCells (vtkUnstructuredGrid *grid, vtkCellArray *elements, vtkCellArray *subelements)
 Create polyhedron in grid from elements and subelements.
 
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkDataSetCreateMesh (const conduit_cpp::Node &topology, const conduit_cpp::Node &coordsets)
 vtkDataSet creation.
 
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkImageDataCreateImageData (const conduit_cpp::Node &coordset)
 Create a vtkImageData from a coordset node.
 
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkRectilinearGridCreateRectilinearGrid (const conduit_cpp::Node &coordset)
 Create a vtkRectilinearGrid from a coordset node.
 
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkStructuredGridCreateStructuredGrid (const conduit_cpp::Node &topology, const conduit_cpp::Node &coordset)
 Create a vtkStructuredGrid from a topology and a coordset nodes.
 
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkDataSetCreateMonoShapedUnstructuredGrid (const conduit_cpp::Node &topologyNode, const conduit_cpp::Node &coordset)
 Create a vtkUnstructuredGrid from a topology and a coordset node.
 
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkDataSetCreateMixedUnstructuredGrid (const conduit_cpp::Node &topologyNode, const conduit_cpp::Node &coords)
 Create a vtkUnstructuredGrid from a coordset and a topology node.
 
VTKIOCATALYSTCONDUIT_EXPORT vtkIdType GetNumberOfPointsInCellType (int vtk_cell_type)
 Return the number of points in VTK cell type.
 
VTKIOCATALYSTCONDUIT_EXPORT int GetCellType (const std::string &shape)
 Get vtk cell type from conduit shape name throw a runtime_error on unsupported type.
 
VTKIOCATALYSTCONDUIT_EXPORT int GetAssociation (const std::string &association)
 Get vtkDataObject attribute type from conduit association string.
 

Function Documentation

◆ FillPartionedDataSet()

VTKIOCATALYSTCONDUIT_EXPORT bool vtkConduitToDataObject::FillPartionedDataSet ( vtkPartitionedDataSet output,
const conduit_cpp::Node &  meshNode 
)

Fill the vtkPartitionedDataSet input.

Create concrete vtkDataSet subclass to set it as partition and add arrays in its DataSetAttributes.

Return true if data was correctly generated, false if an error occurred. Do not throw errors.

◆ FillAMRMesh()

VTKIOCATALYSTCONDUIT_EXPORT bool vtkConduitToDataObject::FillAMRMesh ( vtkOverlappingAMR amr,
const conduit_cpp::Node &  node 
)

Fill the vtkOverlappingAMR input.

◆ CreateMesh()

VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkDataSet > vtkConduitToDataObject::CreateMesh ( const conduit_cpp::Node &  topology,
const conduit_cpp::Node &  coordsets 
)

vtkDataSet creation.

Create a vtkDataSet concrete subclass from the given a topology and a coordsets nodes.

Throw runtime_error on unsupported input.

◆ CreateImageData()

VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkImageData > vtkConduitToDataObject::CreateImageData ( const conduit_cpp::Node &  coordset)

Create a vtkImageData from a coordset node.

◆ CreateRectilinearGrid()

VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkRectilinearGrid > vtkConduitToDataObject::CreateRectilinearGrid ( const conduit_cpp::Node &  coordset)

Create a vtkRectilinearGrid from a coordset node.

◆ CreateStructuredGrid()

VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkStructuredGrid > vtkConduitToDataObject::CreateStructuredGrid ( const conduit_cpp::Node &  topology,
const conduit_cpp::Node &  coordset 
)

Create a vtkStructuredGrid from a topology and a coordset nodes.

◆ CreateMonoShapedUnstructuredGrid()

VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkDataSet > vtkConduitToDataObject::CreateMonoShapedUnstructuredGrid ( const conduit_cpp::Node &  topologyNode,
const conduit_cpp::Node &  coordset 
)

Create a vtkUnstructuredGrid from a topology and a coordset node.

Topology should have a unique cell type, i.e. its "elements/shape" should not be "mixed". see CreateMixedUnstructuredGrid.

◆ CreateMixedUnstructuredGrid()

VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkDataSet > vtkConduitToDataObject::CreateMixedUnstructuredGrid ( const conduit_cpp::Node &  topologyNode,
const conduit_cpp::Node &  coords 
)

Create a vtkUnstructuredGrid from a coordset and a topology node.

Topology "elements/shape" is expected to be "mixed". see CreateMonoShapedUnstructuredGrid throw a runtime_error on invalid node

◆ AddFieldData()

VTKIOCATALYSTCONDUIT_EXPORT bool vtkConduitToDataObject::AddFieldData ( vtkDataObject output,
const conduit_cpp::Node &  stateFields,
bool  isAMReX = false 
)

Add FieldData arrays to output data object.

Return true if node was correctly parsed, false if a fatal error occurred. If isAMReX, data array is added as a CellData

◆ CreatePoints()

VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkPoints > vtkConduitToDataObject::CreatePoints ( const conduit_cpp::Node &  coords)

Create a vtkPoints from a coordset node that respect the following requirements:

  • "type" should be "explicit"
  • "values" should have at max 3 components throw a runtime_error on invalid node

◆ SetPolyhedralCells()

VTKIOCATALYSTCONDUIT_EXPORT void vtkConduitToDataObject::SetPolyhedralCells ( vtkUnstructuredGrid grid,
vtkCellArray elements,
vtkCellArray subelements 
)

Create polyhedron in grid from elements and subelements.

◆ GetNumberOfPointsInCellType()

VTKIOCATALYSTCONDUIT_EXPORT vtkIdType vtkConduitToDataObject::GetNumberOfPointsInCellType ( int  vtk_cell_type)

Return the number of points in VTK cell type.

throw a runtime_error on unsupported type.

◆ GetCellType()

VTKIOCATALYSTCONDUIT_EXPORT int vtkConduitToDataObject::GetCellType ( const std::string &  shape)

Get vtk cell type from conduit shape name throw a runtime_error on unsupported type.

◆ GetAssociation()

VTKIOCATALYSTCONDUIT_EXPORT int vtkConduitToDataObject::GetAssociation ( const std::string &  association)

Get vtkDataObject attribute type from conduit association string.

Supports only "element" and "vertex".

Throw runtime_error on unsupported association.