BioDynaMo  v1.05.120-25dc9790
parallel_vtu_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 // std
17 #include <fstream>
18 #include <sstream>
19 // Paraview
20 #include <vtkXMLPUnstructuredGridWriter.h>
21 #include <vtkXMLUnstructuredGridWriter.h>
22 // BioDynaMo
23 #include "core/param/param.h"
24 #include "core/simulation.h"
25 #include "core/util/string.h"
26 #include "core/util/thread_info.h"
27 
28 namespace bdm {
29 
30 // -----------------------------------------------------------------------------
31 void FixPvtu(const std::string& filename, const std::string& file_prefix,
32  uint64_t pieces) {
33  // read whole pvtu file into buffer
34  std::ifstream ifs(filename);
35  ifs.seekg(0, std::ios::end);
36  size_t size = ifs.tellg();
37  std::string buffer(size, ' ');
38  ifs.seekg(0);
39  ifs.read(&buffer[0], size);
40 
41  // create new file
42  std::string find = Concat("<Piece Source=\"", file_prefix, "_0.vtu\"/>");
43  std::stringstream new_file;
44  auto pos = buffer.find(find);
45  if (pos == std::string::npos) {
46  // No agents in of this type in this timestep
47  return;
48  }
49  new_file << buffer.substr(0, pos);
50  for (uint64_t i = 0; i < pieces; ++i) {
51  new_file << "<Piece Source=\"" << file_prefix << "_" << i << ".vtu\"/>\n";
52  }
53  new_file << buffer.substr(pos + find.size(), buffer.size());
54  // write new file
55  std::ofstream ofs(filename);
56  ofs << new_file.str();
57 }
58 
59 // -----------------------------------------------------------------------------
61  const std::string& folder, const std::string& file_prefix,
62  const std::vector<vtkUnstructuredGrid*>& grids) const {
63  auto* tinfo = ThreadInfo::GetInstance();
64  auto* param = Simulation::GetActive()->GetParam();
65 
66 #pragma omp parallel for schedule(static, 1)
67  for (int i = 0; i < tinfo->GetMaxThreads(); ++i) {
68  if (i == 0) {
69  vtkNew<vtkXMLPUnstructuredGridWriter> pvtu_writer;
70  auto filename = Concat(folder, "/", file_prefix, ".pvtu");
71  pvtu_writer->SetFileName(filename.c_str());
72  auto max_threads = ThreadInfo::GetInstance()->GetMaxThreads();
73  pvtu_writer->SetInputData(grids[0]);
74  pvtu_writer->SetDataModeToBinary();
75  pvtu_writer->SetEncodeAppendedData(false);
76  if (!param->visualization_compress_pv_files) {
77  pvtu_writer->SetCompressorTypeToNone();
78  }
79  pvtu_writer->Write();
80 
81  FixPvtu(filename, file_prefix, max_threads);
82  } else {
83  vtkNew<vtkXMLUnstructuredGridWriter> vtu_writer;
84  vtu_writer->SetGlobalWarningDisplay(false);
85  auto filename = Concat(folder, "/", file_prefix, "_", i, ".vtu");
86  vtu_writer->SetFileName(filename.c_str());
87  vtu_writer->SetInputData(grids[i]);
88  vtu_writer->SetDataModeToBinary();
89  vtu_writer->SetEncodeAppendedData(false);
90  if (!param->visualization_compress_pv_files) {
91  vtu_writer->SetCompressorTypeToNone();
92  }
93  vtu_writer->Write();
94  }
95  }
96 }
97 
98 } // namespace bdm
parallel_vtu_writer.h
bdm::ThreadInfo::GetInstance
static ThreadInfo * GetInstance()
Definition: thread_info.cc:21
bdm
Definition: agent.cc:39
string.h
thread_info.h
bdm::ParallelVtuWriter::operator()
void operator()(const std::string &folder, const std::string &file_prefix, const std::vector< vtkUnstructuredGrid * > &grids) const
Definition: parallel_vtu_writer.cc:60
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
simulation.h
bdm::Simulation::GetParam
const Param * GetParam() const
Returns the simulation parameters.
Definition: simulation.cc:254
bdm::ThreadInfo::GetMaxThreads
int GetMaxThreads() const
Return the maximum number of threads.
Definition: thread_info.h:66
param.h
bdm::Simulation::GetActive
static Simulation * GetActive()
This function returns the currently active Simulation simulation.
Definition: simulation.cc:68
bdm::FixPvtu
void FixPvtu(const std::string &filename, const std::string &file_prefix, uint64_t pieces)
Definition: parallel_vtu_writer.cc:31