OGS
Adaptor.cpp
Go to the documentation of this file.
1 
10 #include "Adaptor.h"
11 
12 #include <vtkCPDataDescription.h>
13 #include <vtkCPInputDataDescription.h>
14 #include <vtkCPProcessor.h>
15 #include <vtkCPPythonScriptPipeline.h>
16 #include <vtkNew.h>
17 #include <vtkUnstructuredGrid.h>
18 
19 #include "MeshLib/Mesh.h"
21 #include "filesystem.h"
22 
23 namespace InSituLib
24 {
25 vtkCPProcessor* Processor = nullptr;
26 
27 void Initialize(BaseLib::ConfigTree const& scripts_config,
28  std::string const& path)
29 {
30  if (Processor == nullptr)
31  {
32  Processor = vtkCPProcessor::New();
33  Processor->Initialize();
34  }
35  else
36  {
37  Processor->RemoveAllPipelines();
38  }
40  for (auto script_config : scripts_config.getConfigSubtreeList("script"))
41  {
42  auto scriptPath =
44  fs::path(script_config.getConfigParameter<std::string>("name"));
45  if (scriptPath.is_relative())
46  {
47  scriptPath = fs::path(path) / scriptPath;
48  }
49  if (!fs::exists(scriptPath))
50  {
51  ERR("In-situ script {:s} does not exist!", scriptPath.string());
52  }
53  INFO("Initializing in-situ script: {:s}", scriptPath.string());
54  vtkNew<vtkCPPythonScriptPipeline> pipeline;
55  pipeline->Initialize(scriptPath.c_str());
56  Processor->AddPipeline(pipeline.GetPointer());
57  }
58 }
59 void Finalize()
60 {
61  if (Processor)
62  {
63  Processor->Delete();
64  Processor = nullptr;
65  }
66 }
67 void CoProcess(MeshLib::Mesh const& mesh, double const time,
68  unsigned int const timeStep, bool const lastTimeStep,
69  std::string output_directory)
70 {
71  if (Processor == nullptr)
72  return;
73 
74  vtkNew<vtkCPDataDescription> dataDescription;
75  dataDescription->AddInput("input");
76  dataDescription->SetTimeData(time, timeStep);
77  if (lastTimeStep == true)
78  {
79  // assume that we want to all the pipelines to execute if it
80  // is the last time step.
81  dataDescription->ForceOutputOn();
82  }
83  if (Processor->RequestDataDescription(dataDescription.GetPointer()) != 0)
84  {
85  INFO("Start InSitu process: timestep #{:d} (t={:g}, last={:d})",
86  timeStep, time, lastTimeStep);
87  vtkNew<MeshLib::VtkMappedMeshSource> vtkSource;
88  vtkSource->SetMesh(&mesh);
89  vtkSource->Update();
90  dataDescription->GetInputDescriptionByName("input")->SetGrid(
91  vtkSource->GetOutput());
92  auto const cwd = fs::current_path();
93  fs::current_path(std::move(output_directory));
94  Processor->CoProcess(dataDescription.GetPointer());
95  fs::current_path(cwd);
96  INFO("End InSitu process.");
97  }
98 }
99 } // namespace InSituLib
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42
Definition of the Mesh class.
VtkMappedMeshSource is a source class to transform OGS meshes into complete vtkUnstructuredGrids....
Range< SubtreeIterator > getConfigSubtreeList(std::string const &root) const
Definition: ConfigTree.cpp:169
void Initialize(BaseLib::ConfigTree const &scripts_config, std::string const &path)
Definition: Adaptor.cpp:27
void Finalize()
Definition: Adaptor.cpp:59
void CoProcess(MeshLib::Mesh const &mesh, double const time, unsigned int const timeStep, bool const lastTimeStep, std::string output_directory)
Definition: Adaptor.cpp:67
vtkCPProcessor * Processor
Definition: Adaptor.cpp:25