10#include "vtkUnstructuredGrid.h"
12#include <vtkCellData.h>
13#include <vtkPointData.h>
14#include <vtkDoubleArray.h>
20 const std::string &file_base_name,
21 const std::vector<std::shared_ptr<const FieldFunctionGridBased>> &ff_list)
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 <<
"\"";
28 throw std::logic_error(fname +
": Cannot be used with empty field-function"
32 const auto& master_ff_ptr = ff_list.front();
33 const auto& master_ff = *master_ff_ptr;
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.");
42 const auto& grid = master_ff.sdm_->Grid();
47 auto cell_data = ugrid->GetCellData();
48 auto point_data = ugrid->GetPointData();
49 for (
const auto& ff_ptr : ff_list)
51 const auto field_vector = ff_ptr->GetGhostedFieldVector();
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();
58 for (
uint c=0; c<num_comps; ++c)
60 std::string component_name = ff_ptr->TextName() +
63 component_name += unknown.component_text_names_[c];
68 point_array->SetName(component_name.c_str());
69 cell_array->SetName(component_name.c_str());
72 for (
const auto& cell : grid.local_cells)
74 const size_t num_nodes = sdm->GetCellNumNodes(cell);
76 if (num_nodes == cell.vertex_ids_.size())
78 double node_average = 0.0;
79 for (
int n=0; n<num_nodes; ++n)
81 const int64_t nmap = sdm->MapDOFLocal(cell,n,uk_man,0,c);
83 const double field_value = field_vector[nmap];
85 point_array->InsertNextValue(field_value);
86 node_average += field_value;
88 node_average /=
static_cast<double>(num_nodes);
89 cell_array->InsertNextValue(node_average);
93 double node_average = 0.0;
94 for (
int n=0; n<num_nodes; ++n)
96 const int64_t nmap = sdm->MapDOFLocal(cell,n,uk_man,0,c);
98 const double field_value = field_vector[nmap];
99 node_average += field_value;
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)
105 point_array->InsertNextValue(node_average);
111 point_data->AddArray(point_array);
112 cell_data->AddArray(cell_array);
118 Chi::log.
Log() <<
"Done exporting field functions to VTK.";
LogStream Log(LOG_LVL level=LOG_0)
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)