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