OGS
InSituLib Namespace Reference

Detailed Description

Functions

void Initialize (BaseLib::ConfigTree const &scripts_config, std::string const &path)
 
void Finalize ()
 
void CoProcess (MeshLib::Mesh const &mesh, double const time, unsigned int const timeStep, bool const lastTimeStep, std::string output_directory)
 

Variables

vtkCPProcessor * Processor = nullptr
 

Function Documentation

◆ CoProcess()

void InSituLib::CoProcess ( MeshLib::Mesh const &  mesh,
double const  time,
unsigned int const  timeStep,
bool const  lastTimeStep,
std::string  output_directory 
)

Definition at line 67 of file Adaptor.cpp.

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 }
void INFO(char const *fmt, Args const &... args)
Definition: Logging.h:32
vtkCPProcessor * Processor
Definition: Adaptor.cpp:25

References INFO(), and Processor.

Referenced by ProcessLib::Output::doOutput(), and ProcessLib::Output::doOutputLastTimestep().

◆ Finalize()

void InSituLib::Finalize ( )

Definition at line 59 of file Adaptor.cpp.

60 {
61  if (Processor)
62  {
63  Processor->Delete();
64  Processor = nullptr;
65  }
66 }

References Processor.

Referenced by main().

◆ Initialize()

void InSituLib::Initialize ( BaseLib::ConfigTree const &  scripts_config,
std::string const &  path 
)
Input File Parameter:
prj__insitu__scripts__script
Input File Parameter:
prj__insitu__scripts__script__name

Definition at line 27 of file Adaptor.cpp.

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 }
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42

References ERR(), BaseLib::ConfigTree::getConfigSubtreeList(), INFO(), and Processor.

Referenced by main(), and MeshLib::VtkMeshNodalCoordinatesTemplate< Scalar >::SetNodes().

Variable Documentation

◆ Processor

vtkCPProcessor* InSituLib::Processor = nullptr

Definition at line 25 of file Adaptor.cpp.

Referenced by CoProcess(), Finalize(), and Initialize().