BioDynaMo  v1.05.120-25dc9790
adaptor.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 
15 #include <cstdlib>
16 #include <fstream>
17 #include <memory>
18 #include <sstream>
19 
24 
25 #ifndef __ROOTCLING__
26 
27 #include <vtkCPDataDescription.h>
28 #include <vtkCPInputDataDescription.h>
29 #include <vtkCPProcessor.h>
30 #include <vtkCPPythonScriptPipeline.h>
31 #include <vtkDoubleArray.h>
32 #include <vtkFieldData.h>
33 #include <vtkIdTypeArray.h>
34 #include <vtkImageData.h>
35 #include <vtkImageDataStreamer.h>
36 #include <vtkIntArray.h>
37 #include <vtkNew.h>
38 #include <vtkPointData.h>
39 #include <vtkPoints.h>
40 #include <vtkStringArray.h>
41 #include <vtkUnstructuredGrid.h>
42 #include <vtkXMLImageDataWriter.h>
43 #include <vtkXMLPImageDataWriter.h>
44 #include <vtkXMLPUnstructuredGridWriter.h>
45 #include <vtkXMLUnstructuredGridWriter.h>
46 
47 namespace bdm {
48 
49 // ----------------------------------------------------------------------------
51  vtkCPProcessor* g_processor_ = nullptr;
52  std::unordered_map<std::string, VtkAgents*> vtk_agents_;
53  std::unordered_map<std::string, VtkDiffusionGrid*> vtk_dgrids_;
54  vtkCPDataDescription* data_description_ = nullptr;
55 };
56 
57 // ----------------------------------------------------------------------------
58 std::atomic<uint64_t> ParaviewAdaptor::counter_;
59 
60 // ----------------------------------------------------------------------------
62  counter_++;
63  impl_ = std::unique_ptr<ParaviewAdaptor::ParaviewImpl>(
65 }
66 
67 // ----------------------------------------------------------------------------
69  auto* param = Simulation::GetActive()->GetParam();
70  counter_--;
71 
72  if (impl_) {
73  if (counter_ == 0 && impl_->g_processor_) {
74  impl_->g_processor_->RemoveAllPipelines();
75  impl_->g_processor_->Finalize();
76  impl_->g_processor_->Delete();
77  impl_->g_processor_ = nullptr;
78  }
79  if (param->export_visualization &&
80  param->visualization_export_generate_pvsm) {
83  }
84 
85  if (impl_->data_description_ != nullptr) {
86  impl_->data_description_->Delete();
87  }
88  for (auto& el : impl_->vtk_agents_) {
89  delete el.second;
90  }
91  for (auto& el : impl_->vtk_dgrids_) {
92  delete el.second;
93  }
94  }
95 }
96 
97 // ----------------------------------------------------------------------------
99  if (!initialized_) {
100  Initialize();
101  initialized_ = true;
102  }
103 
104  auto* sim = Simulation::GetActive();
105  auto* param = sim->GetParam();
106  uint64_t total_steps = sim->GetScheduler()->GetSimulatedSteps();
107  if (total_steps % param->visualization_interval != 0) {
108  return;
109  }
110 
111  real_t time = param->simulation_time_step * total_steps;
112  impl_->data_description_->SetTimeData(time, total_steps);
113 
115 
116  if (param->insitu_visualization) {
118  }
119  if (param->export_visualization) {
121  }
122 }
123 
124 // ----------------------------------------------------------------------------
126  auto* sim = Simulation::GetActive();
127  auto* param = sim->GetParam();
128 
129  if (param->insitu_visualization && impl_->g_processor_ == nullptr) {
130 #ifdef __APPLE__
131  Log::Warning("ParaviewAdaptor",
132  "Insitu visualization is currently not supported on MacOS. "
133  "Please use export visualization.");
134 #endif // __APPLE__
135  impl_->g_processor_ = vtkCPProcessor::New();
136  impl_->g_processor_->Initialize();
137  }
138 
139  if (param->insitu_visualization) {
140  const std::string& script =
141  ParaviewAdaptor::BuildPythonScriptString(param->pv_insitu_pipeline);
142  std::ofstream ofs;
143  auto* sim = Simulation::GetActive();
144  std::string final_python_script_name =
145  Concat(sim->GetOutputDir(), "/insitu_pipeline.py");
146  ofs.open(final_python_script_name);
147  ofs << script;
148  ofs.close();
149  vtkNew<vtkCPPythonScriptPipeline> pipeline;
150  pipeline->Initialize(final_python_script_name.c_str());
151  impl_->g_processor_->AddPipeline(pipeline.GetPointer());
152  }
153 
154  if (impl_->data_description_ == nullptr) {
155  impl_->data_description_ = vtkCPDataDescription::New();
156  } else {
157  impl_->data_description_->Delete();
158  impl_->data_description_ = vtkCPDataDescription::New();
159  }
160  impl_->data_description_->SetTimeData(0, 0);
161 
162  for (auto& pair : param->visualize_agents) {
163  impl_->vtk_agents_[pair.first.c_str()] =
164  new VtkAgents(pair.first.c_str(), impl_->data_description_);
165  }
166  for (auto& entry : param->visualize_diffusion) {
167  impl_->vtk_dgrids_[entry.name] =
168  new VtkDiffusionGrid(entry.name, impl_->data_description_);
169  }
170 }
171 
172 // ----------------------------------------------------------------------------
174 #ifdef BDM_PV_EXPERIMENTAL
175  vtkCommand::kOpenGLCacheEnable = true;
176  vtkCommand::kOpenGLCacheIndex = 0;
177 #endif // BDM_PV_EXPERIMENTAL
178  impl_->g_processor_->RequestDataDescription(impl_->data_description_);
179  impl_->data_description_->ForceOutputOn();
180  impl_->g_processor_->CoProcess(impl_->data_description_);
181 }
182 
183 // ----------------------------------------------------------------------------
186 
187  auto step = impl_->data_description_->GetTimeStep();
188 
189  for (auto& el : impl_->vtk_agents_) {
190  el.second->WriteToFile(step);
191  }
192 
193  for (auto& el : impl_->vtk_dgrids_) {
194  el.second->WriteToFile(step);
195  }
196 }
197 
198 // ----------------------------------------------------------------------------
202  if (impl_->data_description_->GetUserData() == nullptr) {
203  vtkNew<vtkStringArray> json;
204  json->SetName("metadata");
205  json->InsertNextValue(
206  GenerateSimulationInfoJson(impl_->vtk_agents_, impl_->vtk_dgrids_));
207  vtkNew<vtkFieldData> field;
208  field->AddArray(json);
209  impl_->data_description_->SetUserData(field);
210  }
211 }
212 
213 // ----------------------------------------------------------------------------
216  for (auto& pair : impl_->vtk_agents_) {
217  const auto& agents = rm->GetTypeIndex()->GetType(pair.second->GetTClass());
218  pair.second->Update(&agents);
219  }
220 }
221 
222 // ----------------------------------------------------------------------------
225 
226  rm->ForEachDiffusionGrid([&](DiffusionGrid* grid) {
227  auto it = impl_->vtk_dgrids_.find(grid->GetContinuumName());
228  if (it != impl_->vtk_dgrids_.end()) {
229  it->second->Update(grid);
230  }
231  });
232 }
233 
234 // ----------------------------------------------------------------------------
236  // create simulation_info.json
238  std::ofstream ofstr;
239  auto* sim = Simulation::GetActive();
240  ofstr.open(Concat(sim->GetOutputDir(), "/", kSimulationInfoJson));
241  ofstr << GenerateSimulationInfoJson(impl_->vtk_agents_, impl_->vtk_dgrids_);
242  ofstr.close();
244  }
245 }
246 
247 // ----------------------------------------------------------------------------
252  auto* sim = Simulation::GetActive();
253  std::stringstream python_cmd;
254  std::string pv_dir = std::getenv("ParaView_DIR");
255  std::string bdmsys = std::getenv("BDMSYS");
256 
257  python_cmd << pv_dir << "/bin/pvbatch " << bdmsys
258  << "/include/core/visualization/paraview/generate_pv_state.py "
259  << sim->GetOutputDir() << "/" << kSimulationInfoJson;
260  int ret_code = system(python_cmd.str().c_str());
261  if (ret_code) {
262  Log::Fatal("ParaviewAdaptor::GenerateParaviewState",
263  "Error during generation of ParaView state\n", "Command\n",
264  python_cmd.str());
265  }
266 }
267 
268 // ----------------------------------------------------------------------------
270  const std::string& python_script) {
271  std::stringstream script;
272 
273  std::ifstream ifs;
274  ifs.open(python_script, std::ifstream::in);
275  if (!ifs.is_open()) {
276  Log::Fatal("ParaviewAdaptor::BuildPythonScriptString",
277  Concat("Python script (", python_script,
278  ") was not found or could not be opened."));
279  }
280  script << ifs.rdbuf();
281  ifs.close();
282 
283  std::string default_python_script =
284  std::string(std::getenv("BDMSYS")) +
285  std::string(
286  "/include/core/visualization/paraview/default_insitu_pipeline.py");
287 
288  std::ifstream ifs_default;
289  ifs_default.open(default_python_script, std::ifstream::in);
290  if (!ifs_default.is_open()) {
291  Log::Fatal("ParaviewAdaptor::BuildPythonScriptString",
292  Concat("Python script (", default_python_script,
293  ") was not found or could not be opened."));
294  }
295  script << std::endl << ifs_default.rdbuf();
296  ifs_default.close();
297  return script.str();
298 }
299 
300 } // namespace bdm
301 
302 #else
303 
304 namespace bdm {
305 
307 
309 
311 
313 
315 
317 
319  const std::string& python_script) {}
320 } // namespace bdm
321 
322 #endif
bdm::VtkDiffusionGrid
Definition: vtk_diffusion_grid.h:35
bdm::ParaviewAdaptor::ExportVisualization
void ExportVisualization()
Definition: adaptor.cc:184
vtk_agents.h
bdm::ParaviewAdaptor::initialized_
bool initialized_
only needed for insitu visualization
Definition: adaptor.h:57
bdm::ParaviewAdaptor::ParaviewImpl::vtk_dgrids_
std::unordered_map< std::string, VtkDiffusionGrid * > vtk_dgrids_
Definition: adaptor.cc:53
adaptor.h
vtk_diffusion_grid.h
bdm
Definition: agent.cc:39
bdm::ResourceManager::GetTypeIndex
const TypeIndex * GetTypeIndex() const
Definition: resource_manager.h:492
bdm::ParaviewAdaptor::Visualize
void Visualize() override
Visualize one timestep based on the configuration in Param
Definition: adaptor.cc:98
bdm::real_t
double real_t
Definition: real_t.h:21
bdm::ParaviewAdaptor::Initialize
void Initialize()
Definition: adaptor.cc:125
bdm::GenerateSimulationInfoJson
std::string GenerateSimulationInfoJson(const std::unordered_map< std::string, VtkAgents * > &vtk_agents, const std::unordered_map< std::string, VtkDiffusionGrid * > &vtk_dgrids)
Definition: helper.cc:24
bdm::ParaviewAdaptor::counter_
static std::atomic< uint64_t > counter_
Definition: adaptor.h:54
bdm::TypeIndex::GetType
const std::vector< Agent * > & GetType(TClass *tclass) const
Definition: type_index.cc:70
bdm::DiffusionGrid
Definition: diffusion_grid.h:90
bdm::ParaviewAdaptor::ParaviewImpl
Definition: adaptor.cc:50
bdm::Continuum::GetContinuumName
const std::string & GetContinuumName() const
Returns the name of the continuum.
Definition: continuum_interface.h:93
bdm::ParaviewAdaptor::~ParaviewAdaptor
~ParaviewAdaptor() override
Definition: adaptor.cc:68
bdm::Log::Warning
static void Warning(const std::string &location, const Args &... parts)
Prints warning message.
Definition: log.h:67
bdm::ParaviewAdaptor::ParaviewImpl::data_description_
vtkCPDataDescription * data_description_
Definition: adaptor.cc:54
bdm::Simulation::GetResourceManager
ResourceManager * GetResourceManager()
Returns the ResourceManager instance.
Definition: simulation.cc:244
bdm::WriteToFile
void WriteToFile(const std::string &filename, const std::string &content)
Definition: io.cc:83
bdm::ParaviewAdaptor::BuildDiffusionGridVTKStructures
void BuildDiffusionGridVTKStructures()
Create the required vtk objects to visualize diffusion grids.
Definition: adaptor.cc:223
bdm::VtkAgents
Definition: vtk_agents.h:36
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::Log::Fatal
static void Fatal(const std::string &location, const Args &... parts)
Prints fatal error message.
Definition: log.h:115
bdm::ParaviewAdaptor::BuildAgentsVTKStructures
void BuildAgentsVTKStructures()
Create the required vtk objects to visualize agents.
Definition: adaptor.cc:214
bdm::ParaviewAdaptor::ParaviewImpl::g_processor_
vtkCPProcessor * g_processor_
Definition: adaptor.cc:51
bdm::ResourceManager::ForEachDiffusionGrid
void ForEachDiffusionGrid(TFunctor &&f) const
Definition: resource_manager.h:232
bdm::ParaviewAdaptor::GenerateParaviewState
static void GenerateParaviewState()
Definition: adaptor.cc:251
bdm::ParaviewAdaptor::impl_
std::unique_ptr< ParaviewImpl > impl_
Definition: adaptor.h:50
bdm::ParaviewAdaptor::InsituVisualization
void InsituVisualization()
Execute the insitu pipelines that were defined in Initialize
Definition: adaptor.cc:173
bdm::ParaviewAdaptor::BuildPythonScriptString
static std::string BuildPythonScriptString(const std::string &python_script)
Definition: adaptor.cc:269
bdm::ParaviewAdaptor::WriteSimulationInfoJsonFile
void WriteSimulationInfoJsonFile()
Definition: adaptor.cc:235
bdm::Simulation::GetParam
const Param * GetParam() const
Returns the simulation parameters.
Definition: simulation.cc:254
bdm::ParaviewAdaptor::simulation_info_json_generated_
bool simulation_info_json_generated_
Definition: adaptor.h:58
helper.h
bdm::Simulation::GetActive
static Simulation * GetActive()
This function returns the currently active Simulation simulation.
Definition: simulation.cc:68
bdm::ParaviewAdaptor::CreateVtkObjects
void CreateVtkObjects()
Creates the VTK objects that represent the agents in ParaView.
Definition: adaptor.cc:199
bdm::ParaviewAdaptor::ParaviewAdaptor
ParaviewAdaptor()
Definition: adaptor.cc:61
bdm::ParaviewAdaptor::ParaviewImpl::vtk_agents_
std::unordered_map< std::string, VtkAgents * > vtk_agents_
Definition: adaptor.cc:52