Chi-Tech
ff_gridbased_03b_exportvtkwithappend.cc
Go to the documentation of this file.
2
3#include "chi_runtime.h"
4#include "chi_log.h"
5
8
10#include "vtkUnstructuredGrid.h"
11
12#include <vtkCellData.h>
13#include <vtkPointData.h>
14#include <vtkDoubleArray.h>
15
16//###################################################################
17/**Export multiple field functions to VTK.*/
20 const std::string &file_base_name,
21 const std::vector<std::shared_ptr<const FieldFunctionGridBased>> &ff_list)
22{
23 const std::string fname = "chi_physics::FieldFunction::ExportMultipleToVTK";
24 Chi::log.Log() << "Exporting field functions to VTK with file base \""
25 << file_base_name << "\"";
26
27 if (ff_list.empty())
28 throw std::logic_error(fname + ": Cannot be used with empty field-function"
29 " list");
30
31 //============================================= Setup master and check slaves
32 const auto& master_ff_ptr = ff_list.front();
33 const auto& master_ff = *master_ff_ptr;
34
35 for (const auto& ff_ptr : ff_list)
36 if (ff_ptr != master_ff_ptr)
37 if (&ff_ptr->sdm_->Grid() != &master_ff_ptr->sdm_->Grid())
38 throw std::logic_error(fname +
39 ": Cannot be used with field functions based on different grids.");
40
41 //============================================= Get grid
42 const auto& grid = master_ff.sdm_->Grid();
43
44 auto ugrid = chi_mesh::PrepareVtkUnstructuredGrid(grid);
45
46 //============================================= Upload cell/point data
47 auto cell_data = ugrid->GetCellData();
48 auto point_data = ugrid->GetPointData();
49 for (const auto& ff_ptr : ff_list)
50 {
51 const auto field_vector = ff_ptr->GetGhostedFieldVector();
52
53 const auto& uk_man = ff_ptr->GetUnknownManager();
54 const auto& unknown = ff_ptr->Unknown();
55 const auto& sdm = ff_ptr->sdm_;
56 const size_t num_comps = unknown.NumComponents();
57
58 for (uint c=0; c<num_comps; ++c)
59 {
60 std::string component_name = ff_ptr->TextName() +
61 unknown.text_name_;
62 if (num_comps > 1)
63 component_name += unknown.component_text_names_[c];
64
65 vtkNew<vtkDoubleArray> point_array;
66 vtkNew<vtkDoubleArray> cell_array;
67
68 point_array->SetName(component_name.c_str());
69 cell_array->SetName(component_name.c_str());
70
71 //Populate the array here
72 for (const auto& cell : grid.local_cells)
73 {
74 const size_t num_nodes = sdm->GetCellNumNodes(cell);
75
76 if (num_nodes == cell.vertex_ids_.size())
77 {
78 double node_average = 0.0;
79 for (int n=0; n<num_nodes; ++n)
80 {
81 const int64_t nmap = sdm->MapDOFLocal(cell,n,uk_man,0,c);
82
83 const double field_value = field_vector[nmap];
84
85 point_array->InsertNextValue(field_value);
86 node_average += field_value;
87 }//for node
88 node_average /= static_cast<double>(num_nodes);
89 cell_array->InsertNextValue(node_average);
90 }
91 else
92 {
93 double node_average = 0.0;
94 for (int n=0; n<num_nodes; ++n)
95 {
96 const int64_t nmap = sdm->MapDOFLocal(cell,n,uk_man,0,c);
97
98 const double field_value = field_vector[nmap];
99 node_average += field_value;
100 }//for node
101 node_average /= static_cast<double>(num_nodes);
102 cell_array->InsertNextValue(node_average);
103 for (int n=0; n<cell.vertex_ids_.size(); ++n)
104 {
105 point_array->InsertNextValue(node_average);
106 }//for vertex
107 }
108
109 }//for cell
110
111 point_data->AddArray(point_array);
112 cell_data->AddArray(cell_array);
113 }//for component
114 }//for ff_ptr
115
116 chi_mesh::WritePVTUFiles(ugrid, file_base_name);
117
118 Chi::log.Log() << "Done exporting field functions to VTK.";
119}
static chi::ChiLog & log
Definition: chi_runtime.h:81
LogStream Log(LOG_LVL level=LOG_0)
Definition: chi_log.cc:35
static void ExportMultipleToVTK(const std::string &file_base_name, const FFList &ff_list)
vtkNew< vtkUnstructuredGrid > PrepareVtkUnstructuredGrid(const chi_mesh::MeshContinuum &grid, bool discontinuous=true)
void WritePVTUFiles(vtkNew< vtkUnstructuredGrid > &ugrid, const std::string &file_base_name)
#define uint