VTK  9.4.20241221
vtkIOSSUtilities.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
3#ifndef vtkIOSSUtilities_h
4#define vtkIOSSUtilities_h
21#include "vtkDoubleArray.h"
22#include "vtkIOSSReader.h"
23#include "vtkLogger.h"
24#include "vtkObject.h"
25#include "vtkSmartPointer.h"
26#include "vtkTypeInt32Array.h"
27#include "vtkTypeInt64Array.h"
28#include "vtkTypeList.h" // Needed for ArrayList definition
29
30// Ioss includes
31#include <vtk_ioss.h>
32// clang-format off
33#include VTK_IOSS(Ioss_Region.h)
34#include VTK_IOSS(Ioss_Transform.h)
35#include VTK_IOSS(Ioss_StructuredBlock.h)
36#include VTK_IOSS(Ioss_SideSet.h)
37#include VTK_IOSS(Ioss_SideBlock.h)
38// clang-format on
39
40#include <cassert>
41#include <set>
42
43VTK_ABI_NAMESPACE_BEGIN
44class vtkCellArray;
45class vtkDataSet;
46VTK_ABI_NAMESPACE_END
47
48namespace vtkIOSSUtilities
49{
50VTK_ABI_NAMESPACE_BEGIN
51
53{
58};
59
63class Cache
64{
65public:
68
73
79
83 void Clear();
84
85 vtkObject* Find(const Ioss::GroupingEntity* entity, const std::string& cachekey) const;
86 void Insert(const Ioss::GroupingEntity* entity, const std::string& cachekey, vtkObject* array);
87
88private:
89 Cache(const Cache&) = delete;
90 void operator=(const Cache&) = delete;
91
92 class CacheInternals;
93 CacheInternals* Internals;
94};
95
101{
102public:
105
109 std::string GetMessages() const;
110
111private:
112 std::ostringstream Stream;
113 std::ostream* DebugStream;
114 std::ostream* WarningStream;
115};
116
117using EntityNameType = std::pair<vtkTypeUInt64, std::string>;
118
126 vtkTypeList::Create<vtkDoubleArray, vtkTypeInt32Array, vtkTypeInt64Array>>::Result;
127
132std::vector<std::pair<int, double>> GetTime(const Ioss::Region* region);
133
141std::string GetSanitizedBlockName(const Ioss::Region* region, const std::string& name);
142
148template <typename EntityType>
149void GetEntityAndFieldNames(const Ioss::Region* region, const std::vector<EntityType*>& entities,
150 std::set<EntityNameType>& entity_names, std::set<std::string>& field_names)
151{
152 for (const auto& entity : entities)
153 {
154 const int64_t id = entity->property_exists("id") ? entity->get_property("id").get_int() : 0;
155 auto name = vtkIOSSUtilities::GetSanitizedBlockName(region, entity->name());
156 entity_names.insert(EntityNameType{ static_cast<vtkTypeUInt64>(id), name });
157
158 Ioss::NameList attributeNames;
159 entity->field_describe(Ioss::Field::TRANSIENT, &attributeNames);
160 entity->field_describe(Ioss::Field::ATTRIBUTE, &attributeNames);
161 std::copy(
162 attributeNames.begin(), attributeNames.end(), std::inserter(field_names, field_names.end()));
163 }
164}
168template <>
169void GetEntityAndFieldNames<Ioss::SideSet>(const Ioss::Region* region,
170 const std::vector<Ioss::SideSet*>& entities, std::set<EntityNameType>& entity_names,
171 std::set<std::string>& field_names);
172
180
189
196vtkSmartPointer<vtkDataArray> GetData(const Ioss::GroupingEntity* entity,
197 const std::string& fieldname, Ioss::Transform* transform = nullptr, Cache* cache = nullptr,
198 const std::string& cachekey = std::string());
199
208int GetCellType(const Ioss::ElementTopology* topology);
209
217const Ioss::ElementTopology* GetElementTopology(int vtk_cell_type);
218
230 Ioss::GroupingEntity* group_entity, int& vtk_topology_type, Cache* cache = nullptr);
231
238 const Ioss::GroupingEntity* group_entity, Cache* cache = nullptr);
239
246bool IsFieldTransient(const Ioss::GroupingEntity* entity, const std::string& fieldname);
247
249
252std::string GetDisplacementFieldName(Ioss::GroupingEntity* nodeblock);
255
261
266DatabaseFormatType DetectType(const std::string& dbaseName);
267
273DatabaseFormatType GetFormat(const Ioss::GroupingEntity* entity);
274
281std::vector<Ioss::StructuredBlock*> GetMatchingStructuredBlocks(
282 Ioss::Region* region, const std::string& blockname);
283
284VTK_ABI_NAMESPACE_END
285}
286
287#endif
288// VTK-HeaderTest-Exclude: vtkIOSSUtilities.h
object to represent cell connectivity
abstract class to specify dataset behavior
Definition vtkDataSet.h:165
void Clear()
Clears the cache.
void ClearUnused()
Removes all cached entries not accessed since most recent call to ResetAccessCounts.
void ResetAccessCounts()
Call this to clear internal count for hits.
vtkObject * Find(const Ioss::GroupingEntity *entity, const std::string &cachekey) const
void Insert(const Ioss::GroupingEntity *entity, const std::string &cachekey, vtkObject *array)
A helper to instantiate on stack to temporarily redirect non-critical messages emanating from IOSS.
std::string GetMessages() const
Provides access to the accumulated messages.
abstract base class for most VTK objects
Definition vtkObject.h:162
Hold a reference to a vtkObjectBase instance.
internal utilities for vtkIOSSReader
void GetEntityAndFieldNames< Ioss::SideSet >(const Ioss::Region *region, const std::vector< Ioss::SideSet * > &entities, std::set< EntityNameType > &entity_names, std::set< std::string > &field_names)
Specialization for Ioss::SideSet (see paraview/paraview#21231).
void InitializeEnvironmentForIOSS()
Must be called before using any Ioss library functions.
std::vector< std::pair< int, double > > GetTime(const Ioss::Region *region)
Reads time / timestep information from a region.
std::vector< Ioss::StructuredBlock * > GetMatchingStructuredBlocks(Ioss::Region *region, const std::string &blockname)
Returns collection of StructuredBlock's matching the selected blockname.
void GetEntityAndFieldNames(const Ioss::Region *region, const std::vector< EntityType * > &entities, std::set< EntityNameType > &entity_names, std::set< std::string > &field_names)
Populates entitySelection with available entity block (or set) names and populates fieldSelection wit...
DatabaseFormatType GetFormat(const Ioss::GroupingEntity *entity)
Given any GroupingEntity pointer, returns the format that the associated database is in.
std::pair< vtkTypeUInt64, std::string > EntityNameType
bool IsFieldTransient(const Ioss::GroupingEntity *entity, const std::string &fieldname)
Returns true if the field is transient.
vtkSmartPointer< vtkPoints > GetMeshModelCoordinates(const Ioss::GroupingEntity *group_entity, Cache *cache=nullptr)
Read points from the group_entity.
Ioss::EntityType GetIOSSEntityType(vtkIOSSReader::EntityType vtk_type)
For the given vtkIOSSReader::EntityType return the corresponding Ioss::EntityType.
const Ioss::ElementTopology * GetElementTopology(int vtk_cell_type)
Returns an Ioss topology element, if possible, given a VTK cell type.
int GetCellType(const Ioss::ElementTopology *topology)
Returns VTK celltype for a Ioss topology element.
vtkSmartPointer< vtkDataArray > CreateArray(const Ioss::Field &field)
Create an array for the given field.
vtkSmartPointer< vtkCellArray > GetConnectivity(Ioss::GroupingEntity *group_entity, int &vtk_topology_type, Cache *cache=nullptr)
Read connectivity information from the group_entity.
std::string GetSanitizedBlockName(const Ioss::Region *region, const std::string &name)
This is primarily intended for CGNS.
std::string GetDisplacementFieldName(Ioss::GroupingEntity *nodeblock)
Finds a displacement field name.
vtkSmartPointer< vtkDataArray > GetData(const Ioss::GroupingEntity *entity, const std::string &fieldname, Ioss::Transform *transform=nullptr, Cache *cache=nullptr, const std::string &cachekey=std::string())
Returns a VTK array for a given field (fieldname) on the chosen block (or set) entity.
DatabaseFormatType DetectType(const std::string &dbaseName)
Given a filename determines and returns the database type.
Remove all duplicate types from TypeList TList, storing the new list in Result.