VTK  9.1.0
vtkIOSSUtilities.h
Go to the documentation of this file.
1/*=========================================================================
2
3 Program: Visualization Toolkit
4 Module: vtkIOSSUtilities.h
5
6 Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
7 All rights reserved.
8 See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
9
10 This software is distributed WITHOUT ANY WARRANTY; without even
11 the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
12 PURPOSE. See the above copyright notice for more information.
13
14=========================================================================*/
30#ifndef vtkIOSSUtilities_h
31#define vtkIOSSUtilities_h
32
34#include "vtkDoubleArray.h"
35#include "vtkIOSSReader.h"
36#include "vtkLogger.h"
37#include "vtkObject.h"
38#include "vtkSmartPointer.h"
39#include "vtkTypeInt32Array.h"
40#include "vtkTypeInt64Array.h"
41#include "vtkTypeList.h" // Needed for ArrayList definition
42
43// Ioss includes
44#include <vtk_ioss.h>
45// clang-format off
46#include VTK_IOSS(Ioss_Region.h)
47#include VTK_IOSS(Ioss_Transform.h)
48#include VTK_IOSS(Ioss_StructuredBlock.h)
49// clang-format on
50
51#include <cassert>
52#include <set>
53
54class vtkCellArray;
55namespace vtkIOSSUtilities
56{
57
59{
64};
65
69class Cache
70{
71public:
74
79
85
89 void Clear();
90
91 vtkObject* Find(const Ioss::GroupingEntity* entity, const std::string& cachekey) const;
92 void Insert(const Ioss::GroupingEntity* entity, const std::string& cachekey, vtkObject* array);
93
94private:
95 Cache(const Cache&) = delete;
96 void operator=(const Cache&) = delete;
97
98 class CacheInternals;
99 CacheInternals* Internals;
100};
101
102using EntityNameType = std::pair<vtkTypeUInt64, std::string>;
103
111 vtkTypeList::Create<vtkDoubleArray, vtkTypeInt32Array, vtkTypeInt64Array>>::Result;
112
117std::set<double> GetTimeValues(const Ioss::Region* region);
118
126std::string GetSanitizedBlockName(const Ioss::Region* region, const std::string& name);
127
133template <typename EntityType>
134void GetEntityAndFieldNames(const Ioss::Region* region, const std::vector<EntityType*>& entities,
135 std::set<EntityNameType>& entity_names, std::set<std::string>& field_names)
136{
137 for (const auto& entity : entities)
138 {
139 const int64_t id = entity->property_exists("id") ? entity->get_property("id").get_int() : 0;
140 auto name = vtkIOSSUtilities::GetSanitizedBlockName(region, entity->name());
141 entity_names.insert(EntityNameType{ static_cast<vtkTypeUInt64>(id), name });
142
143 Ioss::NameList attributeNames;
144 entity->field_describe(Ioss::Field::TRANSIENT, &attributeNames);
145 entity->field_describe(Ioss::Field::ATTRIBUTE, &attributeNames);
146 std::copy(
147 attributeNames.begin(), attributeNames.end(), std::inserter(field_names, field_names.end()));
148 }
149}
150
158
167
174vtkSmartPointer<vtkDataArray> GetData(const Ioss::GroupingEntity* entity,
175 const std::string& fieldname, Ioss::Transform* transform = nullptr, Cache* cache = nullptr,
176 const std::string& cachekey = std::string());
177
186int GetCellType(const Ioss::ElementTopology* topology);
187
199 Ioss::GroupingEntity* group_entity, int& vtk_topology_type, Cache* cache = nullptr);
200
207 const Ioss::GroupingEntity* group_entity, Cache* cache = nullptr);
208
215bool IsFieldTransient(Ioss::GroupingEntity* entity, const std::string& fieldname);
216
220std::string GetDisplacementFieldName(Ioss::GroupingEntity* nodeblock);
221
227
232DatabaseFormatType DetectType(const std::string& dbaseName);
233
239DatabaseFormatType GetFormat(const Ioss::GroupingEntity* entity);
240
247std::vector<Ioss::StructuredBlock*> GetMatchingStructuredBlocks(
248 Ioss::Region* region, const std::string& blockname);
249};
250
251#endif
252// VTK-HeaderTest-Exclude: vtkIOSSUtilities.h
object to represent cell connectivity
Definition: vtkCellArray.h:290
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)
abstract base class for most VTK objects
Definition: vtkObject.h:82
Hold a reference to a vtkObjectBase instance.
internal utilities for vtkIOSSReader
void InitializeEnvironmentForIOSS()
Must be called before using any Ioss library functions.
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
std::set< double > GetTimeValues(const Ioss::Region *region)
Reads time / timestep information from a region.
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.
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.
bool IsFieldTransient(Ioss::GroupingEntity *entity, const std::string &fieldname)
Returns true if the field is transient.
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.
Definition: vtkTypeList.h:125