VTK  9.4.20250211
vtkConduitToDataObject.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-License-Identifier: BSD-3-Clause
8#ifndef vtkConduitToDataObject_h
9#define vtkConduitToDataObject_h
10
11#include "vtkIOCatalystConduitModule.h" // For windows import/export of shared libraries
12
13#include "vtkDeprecation.h" // for VTK_DEPRECATED_IN_9_5_0
14#include "vtkObject.h" // for ABI namespace
15#include "vtkSmartPointer.h" // for vtkSmartPointer
16
17namespace conduit_cpp
18{
19class Node;
20}
21
22VTK_ABI_NAMESPACE_BEGIN
24class vtkCellArray;
25class vtkDataObject;
26class vtkDataSet;
27class vtkImageData;
30class vtkPoints;
34VTK_ABI_NAMESPACE_END
35
37{
38VTK_ABI_NAMESPACE_BEGIN
39
41
49VTKIOCATALYSTCONDUIT_EXPORT bool FillPartitionedDataSet(
50 vtkPartitionedDataSet* output, const conduit_cpp::Node& meshNode);
51VTK_DEPRECATED_IN_9_5_0("Use FillPartitionedDataSet instead.")
52VTKIOCATALYSTCONDUIT_EXPORT bool FillPartionedDataSet(
53 vtkPartitionedDataSet* output, const conduit_cpp::Node& meshNode);
55
59VTKIOCATALYSTCONDUIT_EXPORT bool FillAMRMesh(vtkOverlappingAMR* amr, const conduit_cpp::Node& node);
60
65
70VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer<vtkDataSet> CreateMesh(
71 const conduit_cpp::Node& topology, const conduit_cpp::Node& coordsets);
72
76VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer<vtkImageData> CreateImageData(
77 const conduit_cpp::Node& coordset);
78
83 const conduit_cpp::Node& coordset);
84
89 const conduit_cpp::Node& topology, const conduit_cpp::Node& coordset);
90
97 const conduit_cpp::Node& topologyNode, const conduit_cpp::Node& coordset);
98
106 const conduit_cpp::Node& topologyNode, const conduit_cpp::Node& coords);
107
109
115VTKIOCATALYSTCONDUIT_EXPORT bool AddFieldData(
116 vtkDataObject* output, const conduit_cpp::Node& stateFields, bool isAMReX = false);
117
124VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer<vtkPoints> CreatePoints(
125 const conduit_cpp::Node& coords);
126
130VTKIOCATALYSTCONDUIT_EXPORT void SetPolyhedralCells(
131 vtkUnstructuredGrid* grid, vtkCellArray* elements, vtkCellArray* subelements);
132
134
138VTKIOCATALYSTCONDUIT_EXPORT vtkIdType GetNumberOfPointsInCellType(int vtk_cell_type);
139
144VTKIOCATALYSTCONDUIT_EXPORT int GetCellType(const std::string& shape);
145
152VTKIOCATALYSTCONDUIT_EXPORT int GetAssociation(const std::string& association);
154
155VTK_ABI_NAMESPACE_END
156}
157
158#endif
159// VTK-HeaderTest-Exclude: vtkConduitToDataObject.h
Abstract superclass for all arrays.
object to represent cell connectivity
general representation of visualization data
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
topologically and geometrically regular array of data
hierarchical dataset of vtkUniformGrids
composite dataset to encapsulates a dataset consisting of partitions.
represent and manipulate 3D points
Definition vtkPoints.h:139
a dataset that is topologically regular with variable spacing in the three coordinate directions
Hold a reference to a vtkObjectBase instance.
topologically regular array of data
dataset represents arbitrary combinations of all possible cell types
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkPoints > CreatePoints(const conduit_cpp::Node &coords)
Create a vtkPoints from a coordset node that respect the following requirements:
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkDataSet > CreateMixedUnstructuredGrid(const conduit_cpp::Node &topologyNode, const conduit_cpp::Node &coords)
Create a vtkUnstructuredGrid from a coordset and a topology node.
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkDataSet > CreateMonoShapedUnstructuredGrid(const conduit_cpp::Node &topologyNode, const conduit_cpp::Node &coordset)
Create a vtkUnstructuredGrid from a topology and a coordset node.
VTKIOCATALYSTCONDUIT_EXPORT bool FillPartionedDataSet(vtkPartitionedDataSet *output, const conduit_cpp::Node &meshNode)
Fill the vtkPartitionedDataSet input.
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkStructuredGrid > CreateStructuredGrid(const conduit_cpp::Node &topology, const conduit_cpp::Node &coordset)
Create a vtkStructuredGrid from a topology and a coordset nodes.
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkImageData > CreateImageData(const conduit_cpp::Node &coordset)
Create a vtkImageData from a coordset node.
VTKIOCATALYSTCONDUIT_EXPORT void SetPolyhedralCells(vtkUnstructuredGrid *grid, vtkCellArray *elements, vtkCellArray *subelements)
Create polyhedron in grid from elements and subelements.
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 bool AddFieldData(vtkDataObject *output, const conduit_cpp::Node &stateFields, bool isAMReX=false)
Add FieldData arrays to output data object.
VTKIOCATALYSTCONDUIT_EXPORT bool FillAMRMesh(vtkOverlappingAMR *amr, const conduit_cpp::Node &node)
Fill the vtkOverlappingAMR input.
VTKIOCATALYSTCONDUIT_EXPORT vtkIdType GetNumberOfPointsInCellType(int vtk_cell_type)
Return the number of points in VTK cell type.
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkRectilinearGrid > CreateRectilinearGrid(const conduit_cpp::Node &coordset)
Create a vtkRectilinearGrid from a coordset node.
VTKIOCATALYSTCONDUIT_EXPORT int GetAssociation(const std::string &association)
Get vtkDataObject attribute type from conduit association string.
VTKIOCATALYSTCONDUIT_EXPORT bool FillPartitionedDataSet(vtkPartitionedDataSet *output, const conduit_cpp::Node &meshNode)
Fill the vtkPartitionedDataSet input.
VTKIOCATALYSTCONDUIT_EXPORT vtkSmartPointer< vtkDataSet > CreateMesh(const conduit_cpp::Node &topology, const conduit_cpp::Node &coordsets)
vtkDataSet creation.
#define VTK_DEPRECATED_IN_9_5_0(reason)
int vtkIdType
Definition vtkType.h:315