VTK  9.5.20250822
vtkIOSSCellGridReaderInternal.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
2// SPDX-FileCopyrightText: Copyright (c) Sandia Corporation
3// SPDX-License-Identifier: BSD-3-Clause
4
6
7#include "vtkCellAttribute.h" // For CellTypeInfo.
8#include "vtkUnstructuredGridFieldAnnotations.h" // For Annotations.
9
10VTK_ABI_NAMESPACE_BEGIN
11
12class vtkDGCell;
14
23{
24public:
26
27 std::vector<vtkSmartPointer<vtkCellGrid>> GetCellGrids(const std::string& blockName,
28 vtkIOSSReader::EntityType vtk_entity_type, const DatabaseHandle& handle, int timestep,
30
31 std::vector<vtkSmartPointer<vtkCellGrid>> GetElementBlock(const std::string& blockName,
32 vtkIOSSReader::EntityType vtk_entity_type, const DatabaseHandle& handle, int timestep,
34
35 std::vector<vtkSmartPointer<vtkCellGrid>> GetElementBlock(const std::string& blockName,
36 const Ioss::GroupingEntity* group_entity, const DatabaseHandle& handle, int timestep,
38
39 std::vector<vtkSmartPointer<vtkCellGrid>> GetSideSet(const std::string& blockName,
40 vtkIOSSReader::EntityType vtk_entity_type, const DatabaseHandle& handle, int timestep,
42
43 std::vector<vtkSmartPointer<vtkCellGrid>> GetNodeSet(const std::string& blockName,
44 vtkIOSSReader::EntityType vtk_entity_type, const DatabaseHandle& handle, int timestep,
46
48 int shape_conn_size, int shape_order, vtkDGCell* dg);
49
51 vtkCellGrid* grid, vtkDGCell* meta, const Ioss::GroupingEntity* group_entity,
52 Ioss::Region* region, const DatabaseHandle& handle, int timestep, bool read_ioss_ids,
53 const std::string& cache_key_suffix);
54
56 vtkCellGrid* grid, vtkDGCell* meta, const Ioss::GroupingEntity* group_entity,
57 Ioss::Region* region, const DatabaseHandle& handle, int timestep, bool read_ioss_ids,
58 const std::string& cache_key_suffix);
59
60 bool ApplyDisplacements(vtkCellGrid* grid, Ioss::Region* region,
61 const Ioss::GroupingEntity* group_entity, const DatabaseHandle& handle, int timestep);
62
64};
65
66VTK_ABI_NAMESPACE_END
Visualization data composed of cells of arbitrary type.
Definition vtkCellGrid.h:49
Base class for a discontinuous Galerkin cells of all shapes.
Definition vtkDGCell.h:44
Store on/off settings for data arrays, etc.
represent and manipulate attribute data in a dataset
Internal methods for the cell-grid version of the IOSS reader.
std::vector< vtkSmartPointer< vtkCellGrid > > GetElementBlock(const std::string &blockName, const Ioss::GroupingEntity *group_entity, const DatabaseHandle &handle, int timestep, vtkIOSSCellGridReader *self)
vtkCellAttribute::CellTypeInfo GetCellGridInfoForBlock(int shape_conn_size, int shape_order, vtkDGCell *dg)
std::vector< vtkSmartPointer< vtkCellGrid > > GetNodeSet(const std::string &blockName, vtkIOSSReader::EntityType vtk_entity_type, const DatabaseHandle &handle, int timestep, vtkIOSSCellGridReader *self)
void GetElementAttributes(vtkDataArraySelection *fieldSelection, vtkDataSetAttributes *arrayGroup, vtkCellGrid *grid, vtkDGCell *meta, const Ioss::GroupingEntity *group_entity, Ioss::Region *region, const DatabaseHandle &handle, int timestep, bool read_ioss_ids, const std::string &cache_key_suffix)
std::vector< vtkSmartPointer< vtkCellGrid > > GetSideSet(const std::string &blockName, vtkIOSSReader::EntityType vtk_entity_type, const DatabaseHandle &handle, int timestep, vtkIOSSCellGridReader *self)
void GetNodalAttributes(vtkDataArraySelection *fieldSelection, vtkDataSetAttributes *arrayGroup, vtkCellGrid *grid, vtkDGCell *meta, const Ioss::GroupingEntity *group_entity, Ioss::Region *region, const DatabaseHandle &handle, int timestep, bool read_ioss_ids, const std::string &cache_key_suffix)
std::vector< vtkSmartPointer< vtkCellGrid > > GetCellGrids(const std::string &blockName, vtkIOSSReader::EntityType vtk_entity_type, const DatabaseHandle &handle, int timestep, vtkIOSSCellGridReader *self)
vtkIOSSCellGridReaderInternal(vtkIOSSCellGridReader *self)
vtkNew< vtkUnstructuredGridFieldAnnotations > Annotations
bool ApplyDisplacements(vtkCellGrid *grid, Ioss::Region *region, const Ioss::GroupingEntity *group_entity, const DatabaseHandle &handle, int timestep)
std::vector< vtkSmartPointer< vtkCellGrid > > GetElementBlock(const std::string &blockName, vtkIOSSReader::EntityType vtk_entity_type, const DatabaseHandle &handle, int timestep, vtkIOSSCellGridReader *self)
Reader for IOSS (Sierra IO System) that produces cell-grid data.
Internal methods and state for the IOSS reader.
Allocate and hold a VTK object.
Definition vtkNew.h:167
std::pair< std::string, int > DatabaseHandle