BioDynaMo  v1.05.120-25dc9790
parallel_vti_writer.cc
Go to the documentation of this file.
1 // -----------------------------------------------------------------------------
2 //
3 // Copyright (C) 2021 CERN & University of Surrey for the benefit of the
4 // BioDynaMo collaboration. All Rights Reserved.
5 //
6 // Licensed under the Apache License, Version 2.0 (the "License");
7 // you may not use this file except in compliance with the License.
8 //
9 // See the LICENSE file distributed with this work for details.
10 // See the NOTICE file distributed with this work for additional information
11 // regarding copyright ownership.
12 //
13 // -----------------------------------------------------------------------------
14 
16 // Paraview
17 #include <vtkDataArray.h>
18 #include <vtkDataCompressor.h>
19 #include <vtkObjectFactory.h>
20 #include <vtkPointData.h>
21 #include <vtkXMLPImageDataWriter.h>
22 // BioDynaMo
23 #include "core/param/param.h"
24 #include "core/simulation.h"
25 #include "core/util/string.h"
26 
27 namespace bdm {
28 
29 // -----------------------------------------------------------------------------
30 VtiWriter::VtiWriter() = default;
31 
32 // -----------------------------------------------------------------------------
33 void VtiWriter::SetWholeExtent(const int* whole_extent) {
34  whole_extent_ = whole_extent;
35 }
36 
37 // -----------------------------------------------------------------------------
39  vtkIndent indent) {
40  this->WriteVectorAttribute("WholeExtent", 6, const_cast<int*>(whole_extent_));
41  vtkImageData* input = this->GetInput();
42  this->WriteVectorAttribute("Origin", 3, input->GetOrigin());
43  this->WriteVectorAttribute("Spacing", 3, input->GetSpacing());
44 }
45 
46 // -----------------------------------------------------------------------------
48 
49 // -----------------------------------------------------------------------------
50 void PvtiWriter::Write(const std::string& folder,
51  const std::string& file_prefix,
52  const std::array<int, 6>& whole_extent,
53  const std::vector<std::array<int, 6>>& piece_extents,
54  vtkImageData* img, VtiWriter* vti) {
55  auto* origin = img->GetOrigin();
56  auto* spacing = img->GetSpacing();
57  auto endianess_str = vti->GetByteOrder() == vtkXMLWriter::LittleEndian
58  ? "LittleEndian"
59  : "BigEndian";
60  auto header_type_str =
61  vti->GetHeaderType() == vtkXMLWriter::UInt32 ? "UInt32" : "UInt64";
62  auto* compressor = vti->GetCompressor();
63  std::string compressor_str =
64  compressor != nullptr ? compressor->GetClassName() : "";
65  // vti->GetDataSetMajorVersion() and vti->GetDataSetMinorVersion() are
66  // protected
67  auto version_str = "0.1";
68 
69  auto filename = Concat(folder, "/", file_prefix, ".pvti");
70  std::ofstream ofs(filename);
71  // header
72  ofs << "<?xml version=\"1.0\"?>\n"
73  << "<VTKFile type=\"PImageData\" "
74  << "version=\"" << version_str << "\" "
75  << "byte_order=\"" << endianess_str << "\" "
76  << "header_type=\"" << header_type_str << "\" "
77  << "compressor=\"" << compressor_str << "\">\n";
78  ofs << "<PImageData WholeExtent=\"" << ArrayToString(whole_extent.data(), 6)
79  << "\" GhostLevel=\"0\" Origin=\"" << ArrayToString(origin, 3)
80  << "\" Spacing=\"" << ArrayToString(spacing, 3) << "\">\n";
81  ofs << " <PPointData>\n";
82 
83  // data arrays
84  auto* pd = img->GetPointData();
85  for (int i = 0; i < pd->GetNumberOfArrays(); ++i) {
86  auto name = pd->GetArray(i)->GetName();
87  auto components = pd->GetArray(i)->GetNumberOfComponents();
88  auto float_size = sizeof(real_t) * 8;
89  ofs << " <PDataArray type=\"Float" << float_size << "\" Name=\""
90  << name << "\" NumberOfComponents=\"" << components << "\"/>\n";
91  }
92  ofs << " </PPointData>\n";
93 
94  // pieces
95  int counter = 0;
96  for (auto& e : piece_extents) {
97  auto vti_filename = Concat(file_prefix, "_", counter++, ".vti");
98  ofs << " <Piece Extent=\"" << ArrayToString(e.data(), 6)
99  << "\" Source=\"" << vti_filename << "\"/>\n";
100  }
101 
102  ofs << " </PImageData>\n";
103  ofs << "</VTKFile>\n";
104 }
105 
106 // -----------------------------------------------------------------------------
108  const std::string& folder, const std::string& file_prefix,
109  const std::vector<vtkImageData*>& images, uint64_t num_pieces,
110  const std::array<int, 6>& whole_extent,
111  const std::vector<std::array<int, 6>>& piece_extents) const {
112  auto* param = Simulation::GetActive()->GetParam();
113 
114 #pragma omp parallel for schedule(static, 1)
115  for (uint64_t i = 0; i < num_pieces; ++i) {
116  auto vti_filename = Concat(folder, "/", file_prefix, "_", i, ".vti");
117  vtkNew<VtiWriter> vti;
118  vti->SetFileName(vti_filename.c_str());
119  vti->SetInputData(images[i]);
120  vti->SetWholeExtent(whole_extent.data());
121  vti->SetDataModeToBinary();
122  vti->SetEncodeAppendedData(false);
123  if (!param->visualization_compress_pv_files) {
124  vti->SetCompressorTypeToNone();
125  }
126  vti->Write();
127 
128  if (i == 0) {
129  PvtiWriter pvti;
130  pvti.Write(folder, file_prefix, whole_extent, piece_extents, images[0],
131  vti);
132  }
133  }
134 }
135 
136 } // namespace bdm
bdm::VtiWriter::whole_extent_
const int * whole_extent_
Definition: parallel_vti_writer.h:44
parallel_vti_writer.h
bdm::VtiWriter
Definition: parallel_vti_writer.h:31
bdm::VtiWriter::SetWholeExtent
void SetWholeExtent(const int *whole_extent)
Definition: parallel_vti_writer.cc:33
bdm
Definition: agent.cc:39
string.h
bdm::PvtiWriter::ArrayToString
std::string ArrayToString(T *data, int length) const
Definition: parallel_vti_writer.h:57
bdm::VtiWriter::VtiWriter
VtiWriter()
bdm::real_t
double real_t
Definition: real_t.h:21
bdm::VtiWriter::WritePrimaryElementAttributes
void WritePrimaryElementAttributes(std::ostream &os, vtkIndent indent) override
Definition: parallel_vti_writer.cc:38
bdm::PvtiWriter
Definition: parallel_vti_writer.h:48
bdm::vtkStandardNewMacro
vtkStandardNewMacro(VtiWriter)
bdm::Concat
std::string Concat(const Args &... parts)
Concatenates all arguments into a string. Equivalent to streaming all arguments into a stringstream a...
Definition: string.h:70
bdm::PvtiWriter::Write
void Write(const std::string &folder, const std::string &file_prefix, const std::array< int, 6 > &whole_extent, const std::vector< std::array< int, 6 >> &piece_extents, vtkImageData *img, VtiWriter *vti)
Definition: parallel_vti_writer.cc:50
bdm::ParallelVtiWriter::operator()
void operator()(const std::string &folder, const std::string &file_prefix, const std::vector< vtkImageData * > &images, uint64_t num_pieces, const std::array< int, 6 > &whole_extent, const std::vector< std::array< int, 6 >> &piece_extents) const
Definition: parallel_vti_writer.cc:107
simulation.h
bdm::Simulation::GetParam
const Param * GetParam() const
Returns the simulation parameters.
Definition: simulation.cc:254
param.h
bdm::Simulation::GetActive
static Simulation * GetActive()
This function returns the currently active Simulation simulation.
Definition: simulation.cc:68