18 #ifndef vtkIOSSUtilities_h
19 #define vtkIOSSUtilities_h
27 #include "vtkTypeInt32Array.h"
28 #include "vtkTypeInt64Array.h"
34 #include VTK_IOSS(Ioss_Region.h)
35 #include VTK_IOSS(Ioss_Transform.h)
36 #include VTK_IOSS(Ioss_StructuredBlock.h)
37 #include VTK_IOSS(Ioss_SideSet.h)
38 #include VTK_IOSS(Ioss_SideBlock.h)
44 VTK_ABI_NAMESPACE_BEGIN
51 VTK_ABI_NAMESPACE_BEGIN
91 void operator=(
const Cache&) =
delete;
94 CacheInternals* Internals;
113 std::ostringstream Stream;
114 std::ostream* DebugStream;
115 std::ostream* WarningStream;
127 vtkTypeList::Create<vtkDoubleArray, vtkTypeInt32Array, vtkTypeInt64Array>>::Result;
133 std::vector<std::pair<int, double>>
GetTime(
const Ioss::Region* region);
149 template <
typename EntityType>
151 std::set<EntityNameType>& entity_names, std::set<std::string>& field_names)
153 for (
const auto& entity : entities)
155 const int64_t
id = entity->property_exists(
"id") ? entity->get_property(
"id").get_int() : 0;
159 Ioss::NameList attributeNames;
160 entity->field_describe(Ioss::Field::TRANSIENT, &attributeNames);
161 entity->field_describe(Ioss::Field::ATTRIBUTE, &attributeNames);
163 attributeNames.begin(), attributeNames.end(), std::inserter(field_names, field_names.end()));
170 void GetEntityAndFieldNames<Ioss::SideSet>(
const Ioss::Region* region,
171 const std::vector<Ioss::SideSet*>& entities, std::set<EntityNameType>& entity_names,
172 std::set<std::string>& field_names);
231 Ioss::GroupingEntity* group_entity,
int& vtk_topology_type,
Cache* cache =
nullptr);
239 const Ioss::GroupingEntity* group_entity,
Cache* cache =
nullptr);
283 Ioss::Region* region,
const std::string& blockname);
285 VTK_ABI_NAMESPACE_END
object to represent cell connectivity
abstract class to specify dataset behavior
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.
void Insert(const Ioss::GroupingEntity *entity, const std::string &cachekey, vtkObject *array)
vtkObject * Find(const Ioss::GroupingEntity *entity, const std::string &cachekey) const
A helper to instantiate on stack to temporarily redirect non-critical messages emanating from IOSS.
CaptureNonErrorMessages()
~CaptureNonErrorMessages()
std::string GetMessages() const
Provides access to the accumulated messages.
abstract base class for most VTK objects
internal utilities for vtkIOSSReader
const Ioss::ElementTopology * GetElementTopology(int vtk_cell_type)
Returns an Ioss topology element, if possible, given a VTK cell type.
vtkSmartPointer< vtkPoints > GetMeshModelCoordinates(const Ioss::GroupingEntity *group_entity, Cache *cache=nullptr)
Read points from the group_entity.
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.
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
Ioss::EntityType GetIOSSEntityType(vtkIOSSReader::EntityType vtk_type)
For the given vtkIOSSReader::EntityType return the corresponding Ioss::EntityType.
std::vector< std::pair< int, double > > GetTime(const Ioss::Region *region)
Reads time / timestep information from a region.
int GetCellType(const Ioss::ElementTopology *topology)
Returns VTK celltype for a Ioss topology element.
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 > 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.
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.