Chi-Tech
chi_grid_vtk_utils_04_writing.cc
Go to the documentation of this file.
2
3#include "chi_meshcontinuum.h"
4
5#include <vtkPoints.h>
6#include <vtkPointData.h>
7#include <vtkCellData.h>
8#include <vtkUnstructuredGrid.h>
9#include <vtkUnsignedIntArray.h>
10
11#include <vtkXMLUnstructuredGridWriter.h>
12#include <vtkXMLPUnstructuredGridWriter.h>
13
14#include "chi_runtime.h"
15#include "chi_log.h"
16
17// ###################################################################
18/**Uploads vertices and cells to an unstructured grid. This routine
19 * also uploads cell material ids (sub-domain ids) and partition ids.*/
22 bool discontinuous /*=true*/)
23{
24 //============================================= Instantiate VTK items
26 vtkNew<vtkPoints> points;
27 vtkNew<vtkIntArray> material_array;
28 vtkNew<vtkUnsignedIntArray> partition_id_array;
29
30 points->SetDataType(VTK_DOUBLE);
31
32 //============================================= Set names
33 material_array->SetName("Material");
34 partition_id_array->SetName("Partition");
35
36 std::vector<uint64_t> vertex_map;
37 if (not discontinuous)
38 {
39 vertex_map.assign(grid.GetGlobalVertexCount(), 0);
40 const size_t num_verts = grid.GetGlobalVertexCount();
41 for (size_t v = 0; v < num_verts; ++v)
42 vertex_map[v] = v;
43 }
44
45 // ############################################# Populate cell information
46 int64_t node_count = 0;
47 for (const auto& cell : grid.local_cells)
48 {
49 if (discontinuous)
51 grid, cell, node_count, points, ugrid);
52 else
53 {
54 for (uint64_t vid : cell.vertex_ids_)
55 {
56 const auto& vertex = grid.vertices[vid];
57 points->InsertNextPoint(vertex.x, vertex.y, vertex.z);
58 vertex_map[vid] = node_count;
59 ++node_count;
60 }
61 chi_mesh::UploadCellGeometryContinuous(cell, vertex_map, ugrid);
62 }
63
64 material_array->InsertNextValue(cell.material_id_);
65 partition_id_array->InsertNextValue(cell.partition_id_);
66 } // for local cells
67 ugrid->SetPoints(points);
68
69 ugrid->GetCellData()->AddArray(material_array);
70 ugrid->GetCellData()->AddArray(partition_id_array);
71
72 return ugrid;
73}
74
75// ###################################################################
76/**Writes an unstructured grid to files (.pvtu and .vtu).*/
78 const std::string& file_base_name)
79{
80 //============================================= Construct file name
81 std::string base_filename = std::string(file_base_name);
82 std::string location_filename = base_filename + std::string("_") +
83 std::to_string(Chi::mpi.location_id) +
84 std::string(".vtu");
85
86 //============================================= Write master file
87 if (Chi::mpi.location_id == 0)
88 {
89 std::string pvtu_file_name = base_filename + std::string(".pvtu");
90
92
93 pgrid_writer->EncodeAppendedDataOff();
94 pgrid_writer->SetFileName(pvtu_file_name.c_str());
95 pgrid_writer->SetNumberOfPieces(Chi::mpi.process_count);
96 pgrid_writer->SetStartPiece(Chi::mpi.location_id);
97 pgrid_writer->SetEndPiece(Chi::mpi.process_count - 1);
98 pgrid_writer->SetInputData(ugrid);
99
100 pgrid_writer->Write();
101 }
103
104 //============================================= Serial output each piece
106
107 grid_writer->SetInputData(ugrid);
108 grid_writer->SetFileName(location_filename.c_str());
109
110 grid_writer->Write();
111}
static chi::MPI_Info & mpi
Definition: chi_runtime.h:78
void Barrier() const
Definition: mpi_info.cc:38
LocalCellHandler local_cells
uint64_t GetGlobalVertexCount() const
vtkNew< vtkUnstructuredGrid > PrepareVtkUnstructuredGrid(const chi_mesh::MeshContinuum &grid, bool discontinuous=true)
void UploadCellGeometryContinuous(const chi_mesh::Cell &cell, const std::vector< uint64_t > &vertex_map, vtkNew< vtkUnstructuredGrid > &ugrid)
void WritePVTUFiles(vtkNew< vtkUnstructuredGrid > &ugrid, const std::string &file_base_name)
void UploadCellGeometryDiscontinuous(const chi_mesh::MeshContinuum &grid, const chi_mesh::Cell &cell, int64_t &node_counter, vtkNew< vtkPoints > &points, vtkNew< vtkUnstructuredGrid > &ugrid)