Chi-Tech
unpartmesh_01c_readfromensight.cc
Go to the documentation of this file.
2
4
5#include "utils/chi_utils.h"
6
7#include <vtkSmartPointer.h>
8#include <vtkUnstructuredGrid.h>
9#include <vtkEnSightGoldBinaryReader.h>
10#include <vtkMultiBlockDataSet.h>
11
12#include <vtkInformation.h>
13
14#include "chi_runtime.h"
15#include "chi_log.h"
16
17#define ErrorReadingFile(fname) \
18std::runtime_error("Failed to open file: " + options.file_name + \
19" in call to " + #fname + ".")
20
21//###################################################################
22/**Reads an Ensight-Gold unstructured mesh.*/
25{
26 Chi::log.Log() << "Reading Ensight-Gold file: " << options.file_name << ".";
27
28 //======================================== Attempt to open file
29 std::ifstream file;
30 file.open(options.file_name);
31 if (!file.is_open()) throw ErrorReadingFile(ReadFromEnsightGold);
32 file.close();
33
34 //======================================== Read the file
35 mesh_options_ = options;
37 reader->SetCaseFileName(options.file_name.c_str());
38
39 if (not reader->CanReadFile(options.file_name.c_str()))
40 throw std::logic_error("Unable to read file-type with this routine");
41 reader->UpdateInformation();
42 reader->Update();
43
44 //======================================== Get all the grid blocks
45 auto multiblock = reader->GetOutput();
46
47 std::vector<vtkUGridPtrAndName> grid_blocks;
48 auto iter_a = multiblock->NewIterator();
49 iter_a->GoToFirstItem();
50 while (not iter_a->IsDoneWithTraversal())
51 {
52 auto block_a = iter_a->GetCurrentDataObject();
53
54 const std::string block_name = chi::StringTrim(
55 iter_a->GetCurrentMetaData()->Get(vtkCompositeDataSet::NAME()));
56
57 if (block_a->GetDataObjectType() == VTK_UNSTRUCTURED_GRID)
58 {
59 grid_blocks.emplace_back(
60 vtkUnstructuredGrid::SafeDownCast(block_a),
62 iter_a->GetCurrentMetaData()->Get(vtkCompositeDataSet::NAME())));
63
65 << "Reading block " << block_name
66 << " Number of cells: " << grid_blocks.back().first->GetNumberOfCells()
67 << " Number of points: " << grid_blocks.back().first->GetNumberOfPoints();
68 }
69
70 iter_a->GoToNextItem();
71 }
72
73 //======================================== Get the main + bndry blocks
74 const int max_dimension = chi_mesh::FindHighestDimension(grid_blocks);
75 Chi::log.Log0Verbose1() << "Maximum dimension : " << max_dimension << "\n";
76 std::vector<vtkUGridPtrAndName> domain_grid_blocks =
77 chi_mesh::GetBlocksOfDesiredDimension(grid_blocks, max_dimension);
78 std::vector<vtkUGridPtrAndName> bndry_grid_blocks =
79 chi_mesh::GetBlocksOfDesiredDimension(grid_blocks, max_dimension-1);
80
81
82
83 //======================================== Process blocks
84 chi_mesh::SetBlockIDArrays(domain_grid_blocks);
85 auto ugrid = chi_mesh::ConsolidateGridBlocks(domain_grid_blocks);
86
87 //======================================== Copy Data
88 // Material-IDs will get set form block-id arrays
89 CopyUGridCellsAndPoints(*ugrid, options.scale, max_dimension);
90
91 //======================================== Always do this
93 switch (max_dimension)
94 {
95 case 1: dimension = DIMENSION_1; break;
96 case 2: dimension = DIMENSION_2; break;
97 case 3: dimension = DIMENSION_3; break;
98 default: break;
99 }
100
101 attributes_ = dimension | UNSTRUCTURED;
102
105
106 //======================================== Set boundary ids
107 SetBoundaryIDsFromBlocks(bndry_grid_blocks);
108
109 Chi::log.Log() << "Done reading Ensight-Gold file: "
110 << options.file_name << ".";
111}
static chi::ChiLog & log
Definition: chi_runtime.h:81
LogStream Log(LOG_LVL level=LOG_0)
Definition: chi_log.cc:35
LogStream Log0Verbose1()
Definition: chi_log.h:234
void CopyUGridCellsAndPoints(vtkUnstructuredGrid &ugrid, double scale, int dimension_to_copy)
void ReadFromEnsightGold(const Options &options)
void SetBoundaryIDsFromBlocks(std::vector< vtkUGridPtrAndName > &bndry_grid_blocks)
vtkUGridPtr ConsolidateGridBlocks(std::vector< vtkUGridPtrAndName > &ugrid_blocks, const std::string &block_id_array_name="BlockID")
MeshAttributes
Definition: chi_mesh.h:70
@ UNSTRUCTURED
Definition: chi_mesh.h:77
@ DIMENSION_1
Definition: chi_mesh.h:72
@ DIMENSION_2
Definition: chi_mesh.h:73
@ NONE
Definition: chi_mesh.h:71
@ DIMENSION_3
Definition: chi_mesh.h:74
void SetBlockIDArrays(std::vector< vtkUGridPtrAndName > &ugrid_blocks)
int FindHighestDimension(std::vector< vtkUGridPtrAndName > &ugrid_blocks)
std::vector< vtkUGridPtrAndName > GetBlocksOfDesiredDimension(std::vector< vtkUGridPtrAndName > &ugrid_blocks, int desired_dimension)
std::string StringTrim(const std::string &s)
#define ErrorReadingFile(fname)